int main(int argc, char *argv[]){

  UserInput_t XLAL_INIT_DECL(uvar);
  static ConfigVariables config;

  /* sft related variables */
  MultiSFTVector *inputSFTs = NULL;
  MultiPSDVector *multiPSDs = NULL;
  MultiNoiseWeights *multiWeights = NULL;
  MultiLIGOTimeGPSVector *multiTimes = NULL;
  MultiLALDetector multiDetectors;
  MultiDetectorStateSeries *multiStates = NULL;
  MultiAMCoeffs *multiCoeffs = NULL;
  SFTIndexList *sftIndices = NULL;
  SFTPairIndexList *sftPairs = NULL;
  REAL8Vector *shiftedFreqs = NULL;
  UINT4Vector *lowestBins = NULL;
  COMPLEX8Vector *expSignalPhases = NULL;
  REAL8VectorSequence *sincList = NULL;
  PulsarDopplerParams XLAL_INIT_DECL(dopplerpos);
  PulsarDopplerParams thisBinaryTemplate, binaryTemplateSpacings;
  PulsarDopplerParams minBinaryTemplate, maxBinaryTemplate;
  SkyPosition XLAL_INIT_DECL(skyPos);
  MultiSSBtimes *multiBinaryTimes = NULL;

  INT4  k;
  UINT4 j;
  REAL8 fMin, fMax; /* min and max frequencies read from SFTs */
  REAL8 deltaF; /* frequency resolution associated with time baseline of SFTs */

  REAL8 diagff = 0; /*diagonal metric components*/
  REAL8 diagaa = 0;
  REAL8 diagTT = 0;
  REAL8 diagpp = 1;
  REAL8 ccStat = 0;
  REAL8 evSquared=0;
  REAL8 estSens=0; /*estimated sensitivity(4.13)*/
  BOOLEAN dopplerShiftFlag = TRUE;
  toplist_t *ccToplist=NULL;
  CrossCorrBinaryOutputEntry thisCandidate;
  UINT4 checksum;

  LogPrintf (LOG_CRITICAL, "Starting time\n"); /*for debug convenience to record calculating time*/
  /* initialize and register user variables */
  LIGOTimeGPS computingStartGPSTime, computingEndGPSTime;
  XLALGPSTimeNow (&computingStartGPSTime); /* record the rough starting GPS time*/

  if ( XLALInitUserVars( &uvar ) != XLAL_SUCCESS ) {
    LogPrintf ( LOG_CRITICAL, "%s: XLALInitUserVars() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* read user input from the command line or config file */
  if ( XLALUserVarReadAllInput ( argc, argv ) != XLAL_SUCCESS ) {
    LogPrintf ( LOG_CRITICAL, "%s: XLALUserVarReadAllInput() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  if (uvar.help)	/* if help was requested, then exit */
    return 0;

  CHAR *VCSInfoString = XLALGetVersionString(0);     /**<LAL + LALapps Vsersion string*/
  /*If the version information was requested, output it and exit*/
  if ( uvar.version ){
    XLAL_CHECK ( VCSInfoString != NULL, XLAL_EFUNC, "XLALGetVersionString(0) failed.\n" );
    printf ("%s\n", VCSInfoString );
    exit (0);
  }

  /* configure useful variables based on user input */
  if ( XLALInitializeConfigVars ( &config, &uvar) != XLAL_SUCCESS ) {
    LogPrintf ( LOG_CRITICAL, "%s: XLALInitUserVars() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  deltaF = config.catalog->data[0].header.deltaF;
  REAL8 Tsft = 1.0 / deltaF;

  if (XLALUserVarWasSet(&uvar.spacingF) && XLALUserVarWasSet(&uvar.mismatchF))
    LogPrintf (LOG_CRITICAL, "spacingF and mismatchF are both set, use spacingF %.9g by default\n\n", uvar.spacingF);
  if (XLALUserVarWasSet(&uvar.spacingA) && XLALUserVarWasSet(&uvar.mismatchA))
    LogPrintf (LOG_CRITICAL, "spacingA and mismatchA are both set, use spacingA %.9g by default\n\n", uvar.spacingA);
  if (XLALUserVarWasSet(&uvar.spacingT) && XLALUserVarWasSet(&uvar.mismatchT))
    LogPrintf (LOG_CRITICAL, "spacingT and mismatchT are both set, use spacingT %.9g by default\n\n", uvar.spacingT);
  if (XLALUserVarWasSet(&uvar.spacingP) && XLALUserVarWasSet(&uvar.mismatchP))
    LogPrintf (LOG_CRITICAL, "spacingP and mismatchP are both set, use spacingP %.9g by default\n\n", uvar.spacingP);

  /* create the toplist */
  create_crossCorrBinary_toplist( &ccToplist, uvar.numCand);
  /* now read the data */

  /* /\* get SFT parameters so that we can initialise search frequency resolutions *\/ */
  /* /\* calculate deltaF_SFT *\/ */
  /* deltaF_SFT = catalog->data[0].header.deltaF;  /\* frequency resolution *\/ */
  /* timeBase= 1.0/deltaF_SFT; /\* sft baseline *\/ */

  /* /\* catalog is ordered in time so we can get start, end time and tObs *\/ */
  /* firstTimeStamp = catalog->data[0].header.epoch; */
  /* lastTimeStamp = catalog->data[catalog->length - 1].header.epoch; */
  /* tObs = XLALGPSDiff( &lastTimeStamp, &firstTimeStamp ) + timeBase; */

  /* /\*set pulsar reference time *\/ */
  /* if (LALUserVarWasSet ( &uvar_refTime )) { */
  /*   XLALGPSSetREAL8(&refTime, uvar_refTime); */
  /* }  */
  /* else {	/\*if refTime is not set, set it to midpoint of sfts*\/ */
  /*   XLALGPSSetREAL8(&refTime, (0.5*tObs) + XLALGPSGetREAL8(&firstTimeStamp));  */
  /* } */

  /* /\* set frequency resolution defaults if not set by user *\/ */
  /* if (!(LALUserVarWasSet (&uvar_fResolution))) { */
  /*   uvar_fResolution = 1/tObs; */
  /* } */

  /* { */
  /*   /\* block for calculating frequency range to read from SFTs *\/ */
  /*   /\* user specifies freq and fdot range at reftime */
  /*      we translate this range of fdots to start and endtime and find */
  /*      the largest frequency band required to cover the  */
  /*      frequency evolution  *\/ */
  /*   PulsarSpinRange spinRange_startTime; /\**< freq and fdot range at start-time of observation *\/ */
  /*   PulsarSpinRange spinRange_endTime;   /\**< freq and fdot range at end-time of observation *\/ */
  /*   PulsarSpinRange spinRange_refTime;   /\**< freq and fdot range at the reference time *\/ */

  /*   REAL8 startTime_freqLo, startTime_freqHi, endTime_freqLo, endTime_freqHi, freqLo, freqHi; */

  /*   REAL8Vector *fdotsMin=NULL; */
  /*   REAL8Vector *fdotsMax=NULL; */

  /*   UINT4 k; */

  /*   fdotsMin = (REAL8Vector *)LALCalloc(1, sizeof(REAL8Vector)); */
  /*   fdotsMin->length = N_SPINDOWN_DERIVS; */
  /*   fdotsMin->data = (REAL8 *)LALCalloc(fdotsMin->length, sizeof(REAL8)); */

  /*   fdotsMax = (REAL8Vector *)LALCalloc(1, sizeof(REAL8Vector)); */
  /*   fdotsMax->length = N_SPINDOWN_DERIVS; */
  /*   fdotsMax->data = (REAL8 *)LALCalloc(fdotsMax->length, sizeof(REAL8)); */

  /*   XLAL_INIT_MEM(spinRange_startTime); */
  /*   XLAL_INIT_MEM(spinRange_endTime); */
  /*   XLAL_INIT_MEM(spinRange_refTime); */

  /*   spinRange_refTime.refTime = refTime; */
  /*   spinRange_refTime.fkdot[0] = uvar_f0; */
  /*   spinRange_refTime.fkdotBand[0] = uvar_fBand; */
  /* } */

  /* FIXME: need to correct fMin and fMax for Doppler shift, rngmedian bins and spindown range */
  /* this is essentially just a place holder for now */
  /* FIXME: this running median buffer is overkill, since the running median block need not be centered on the search frequency */
  REAL8 vMax = LAL_TWOPI * (uvar.orbitAsiniSec + uvar.orbitAsiniSecBand) / uvar.orbitPSec + LAL_TWOPI * LAL_REARTH_SI / (LAL_DAYSID_SI * LAL_C_SI) + LAL_TWOPI * LAL_AU_SI/(LAL_YRSID_SI * LAL_C_SI); /*calculate the maximum relative velocity in speed of light*/
  fMin = uvar.fStart * (1 - vMax) - 0.5 * uvar.rngMedBlock * deltaF;
  fMax = (uvar.fStart + uvar.fBand) * (1 + vMax) + 0.5 * uvar.rngMedBlock * deltaF;

  /* read the SFTs*/
  if ((inputSFTs = XLALLoadMultiSFTs ( config.catalog, fMin, fMax)) == NULL){
    LogPrintf ( LOG_CRITICAL, "%s: XLALLoadMultiSFTs() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* calculate the psd and normalize the SFTs */
  if (( multiPSDs =  XLALNormalizeMultiSFTVect ( inputSFTs, uvar.rngMedBlock, NULL )) == NULL){
    LogPrintf ( LOG_CRITICAL, "%s: XLALNormalizeMultiSFTVect() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* compute the noise weights for the AM coefficients */
  if (( multiWeights = XLALComputeMultiNoiseWeights ( multiPSDs, uvar.rngMedBlock, 0 )) == NULL){
    LogPrintf ( LOG_CRITICAL, "%s: XLALComputeMultiNoiseWeights() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* read the timestamps from the SFTs */
  if ((multiTimes = XLALExtractMultiTimestampsFromSFTs ( inputSFTs )) == NULL){
    LogPrintf ( LOG_CRITICAL, "%s: XLALExtractMultiTimestampsFromSFTs() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* read the detector information from the SFTs */
  if ( XLALMultiLALDetectorFromMultiSFTs ( &multiDetectors, inputSFTs ) != XLAL_SUCCESS){
    LogPrintf ( LOG_CRITICAL, "%s: XLALMultiLALDetectorFromMultiSFTs() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* Find the detector state for each SFT */
  /* Offset by Tsft/2 to get midpoint as timestamp */
  if ((multiStates = XLALGetMultiDetectorStates ( multiTimes, &multiDetectors, config.edat, 0.5 * Tsft )) == NULL){
    LogPrintf ( LOG_CRITICAL, "%s: XLALGetMultiDetectorStates() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* Note this is specialized to a single sky position */
  /* This might need to be moved into the config variables */
  skyPos.system = COORDINATESYSTEM_EQUATORIAL;
  skyPos.longitude = uvar.alphaRad;
  skyPos.latitude  = uvar.deltaRad;

  /* Calculate the AM coefficients (a,b) for each SFT */
  if ((multiCoeffs = XLALComputeMultiAMCoeffs ( multiStates, multiWeights, skyPos )) == NULL){
    LogPrintf ( LOG_CRITICAL, "%s: XLALComputeMultiAMCoeffs() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* Construct the flat list of SFTs (this sort of replicates the
     catalog, but there's not an obvious way to get the information
     back) */

  if ( ( XLALCreateSFTIndexListFromMultiSFTVect( &sftIndices, inputSFTs ) != XLAL_SUCCESS ) ) {
    LogPrintf ( LOG_CRITICAL, "%s: XLALCreateSFTIndexListFromMultiSFTVect() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* Construct the list of SFT pairs */
#define PCC_SFTPAIR_HEADER "# The length of SFT-pair list is %u #\n"
#define PCC_SFTPAIR_BODY "%u %u\n"
#define PCC_SFT_HEADER "# The length of SFT list is %u #\n"
#define PCC_SFT_BODY "%s %d %d\n"
  FILE *fp = NULL;

  if (XLALUserVarWasSet(&uvar.pairListInputFilename)) { /* If the user provided a list for reading, use it */
    if((sftPairs = XLALCalloc(1, sizeof(sftPairs))) == NULL){
      XLAL_ERROR(XLAL_ENOMEM);
    }
    if((fp = fopen(uvar.pairListInputFilename, "r")) == NULL){
      LogPrintf ( LOG_CRITICAL, "didn't find SFT-pair list file with given input name\n");
      XLAL_ERROR( XLAL_EFUNC );
    }
    if(fscanf(fp,PCC_SFTPAIR_HEADER,&sftPairs->length)==EOF){
      LogPrintf ( LOG_CRITICAL, "can't read the length of SFT-pair list from the header\n");
      XLAL_ERROR( XLAL_EFUNC );
    }

    if((sftPairs->data = XLALCalloc(sftPairs->length, sizeof(*sftPairs->data)))==NULL){
      XLALFree(sftPairs);
      XLAL_ERROR(XLAL_ENOMEM);
    }

    for(j = 0; j < sftPairs->length; j++){ /*read in  the SFT-pair list */
      if(fscanf(fp,PCC_SFTPAIR_BODY, &sftPairs->data[j].sftNum[0], &sftPairs->data[j].sftNum[1])==EOF){
	LogPrintf ( LOG_CRITICAL, "The length of SFT-pair list doesn't match!");
	XLAL_ERROR( XLAL_EFUNC );
      }
    }
    fclose(fp);

  }

  else { /* if not, construct the list of pairs */
    if ( ( XLALCreateSFTPairIndexList( &sftPairs, sftIndices, inputSFTs, uvar.maxLag, uvar.inclAutoCorr ) != XLAL_SUCCESS ) ) {
      LogPrintf ( LOG_CRITICAL, "%s: XLALCreateSFTPairIndexList() failed with errno=%d\n", __func__, xlalErrno );
      XLAL_ERROR( XLAL_EFUNC );
    }
  }

  if (XLALUserVarWasSet(&uvar.pairListOutputFilename)) { /* Write the list of pairs to a file, if a name was provided */
    if((fp = fopen(uvar.pairListOutputFilename, "w")) == NULL){
      LogPrintf ( LOG_CRITICAL, "Can't write in SFT-pair list \n");
      XLAL_ERROR( XLAL_EFUNC );
    }
    fprintf(fp,PCC_SFTPAIR_HEADER, sftPairs->length ); /*output the length of SFT-pair list to the header*/
    for(j = 0; j < sftPairs->length; j++){
      fprintf(fp,PCC_SFTPAIR_BODY, sftPairs->data[j].sftNum[0], sftPairs->data[j].sftNum[1]);
    }
    fclose(fp);
  }

  if (XLALUserVarWasSet(&uvar.sftListOutputFilename)) { /* Write the list of SFTs to a file for sanity-checking purposes */
    if((fp = fopen(uvar.sftListOutputFilename, "w")) == NULL){
      LogPrintf ( LOG_CRITICAL, "Can't write in flat SFT list \n");
      XLAL_ERROR( XLAL_EFUNC );
    }
    fprintf(fp,PCC_SFT_HEADER, sftIndices->length ); /*output the length of SFT list to the header*/
    for(j = 0; j < sftIndices->length; j++){ /*output the SFT list */
      fprintf(fp,PCC_SFT_BODY, inputSFTs->data[sftIndices->data[j].detInd]->data[sftIndices->data[j].sftInd].name, inputSFTs->data[sftIndices->data[j].detInd]->data[sftIndices->data[j].sftInd].epoch.gpsSeconds, inputSFTs->data[sftIndices->data[j].detInd]->data[sftIndices->data[j].sftInd].epoch.gpsNanoSeconds);
    }
    fclose(fp);
  }

  else if(XLALUserVarWasSet(&uvar.sftListInputFilename)){ /*do a sanity check of the order of SFTs list if the name of input SFT list is given*/
    UINT4 numofsft=0;
    if((fp = fopen(uvar.sftListInputFilename, "r")) == NULL){
      LogPrintf ( LOG_CRITICAL, "Can't read in flat SFT list \n");
      XLAL_ERROR( XLAL_EFUNC );
    }
    if (fscanf(fp, PCC_SFT_HEADER, &numofsft)==EOF){
      LogPrintf ( LOG_CRITICAL, "can't read in the length of SFT list from header\n");
      XLAL_ERROR( XLAL_EFUNC );
    }

    CHARVectorSequence *checkDet=NULL;
    if ((checkDet = XLALCreateCHARVectorSequence (numofsft, LALNameLength) ) == NULL){
      LogPrintf ( LOG_CRITICAL, "%s: XLALCreateCHARVector() failed with errno=%d\n", __func__, xlalErrno );
      XLAL_ERROR( XLAL_EFUNC );
    }
    INT4 checkGPS[numofsft], checkGPSns[numofsft];
    if(numofsft == sftIndices->length){
      for (j=0; j<numofsft; j++){
	if( fscanf(fp,PCC_SFT_BODY,&checkDet->data[j * LALNameLength], &checkGPS[j], &checkGPSns[j])==EOF){
	  LogPrintf ( LOG_CRITICAL, "The length of SFT list doesn't match\n");
	  XLAL_ERROR( XLAL_EFUNC );
	}
	if(strcmp( inputSFTs->data[sftIndices->data[j].detInd]->data[sftIndices->data[j].sftInd].name, &checkDet->data[j * LALNameLength] ) != 0
	   ||inputSFTs->data[sftIndices->data[j].detInd]->data[sftIndices->data[j].sftInd].epoch.gpsSeconds != checkGPS[j]
	   ||inputSFTs->data[sftIndices->data[j].detInd]->data[sftIndices->data[j].sftInd].epoch.gpsNanoSeconds != checkGPSns[j] ){
	  LogPrintf ( LOG_CRITICAL, "The order of SFTs has been changed, it's the end of civilization\n");
	  XLAL_ERROR( XLAL_EFUNC );
	}
      }
      fclose(fp);
      XLALDestroyCHARVectorSequence(checkDet);
    }
    else{
      LogPrintf ( LOG_CRITICAL, "Run for your life, the length of SFT list doesn't match");
      XLAL_ERROR( XLAL_EFUNC );
    }
  }
  else
    {

    }

  /* Get weighting factors for calculation of metric */
  /* note that the sigma-squared is now absorbed into the curly G
     because the AM coefficients are noise-weighted. */
  REAL8Vector *GammaAve = NULL;
  REAL8Vector *GammaCirc = NULL;
  if ( ( XLALCalculateCrossCorrGammas( &GammaAve, &GammaCirc, sftPairs, sftIndices, multiCoeffs)  != XLAL_SUCCESS ) ) {
    LogPrintf ( LOG_CRITICAL, "%s: XLALCalculateCrossCorrGammas() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

#define PCC_GAMMA_HEADER "# The normalization Sinv_Tsft is %g #\n"
#define PCC_GAMMA_BODY "%.10g\n"
  if (XLALUserVarWasSet(&uvar.gammaAveOutputFilename)) { /* Write the aa+bb weight for each pair to a file, if a name was provided */
    if((fp = fopen(uvar.gammaAveOutputFilename, "w")) == NULL) {
      LogPrintf ( LOG_CRITICAL, "Can't write in Gamma_ave list \n");
      XLAL_ERROR( XLAL_EFUNC );
    }
    fprintf(fp,PCC_GAMMA_HEADER, multiWeights->Sinv_Tsft); /*output the normalization factor to the header*/
    for(j = 0; j < sftPairs->length; j++){
      fprintf(fp,PCC_GAMMA_BODY, GammaAve->data[j]);
    }
    fclose(fp);
  }
  if (XLALUserVarWasSet(&uvar.gammaCircOutputFilename)) { /* Write the ab-ba weight for each pair to a file, if a name was provided */
    if((fp = fopen(uvar.gammaCircOutputFilename, "w")) == NULL) {
      LogPrintf ( LOG_CRITICAL, "Can't write in Gamma_circ list \n");
      XLAL_ERROR( XLAL_EFUNC );
    }
    fprintf(fp,PCC_GAMMA_HEADER, multiWeights->Sinv_Tsft); /*output the normalization factor to the header*/
    for(j = 0; j < sftPairs->length; j++){
      fprintf(fp,PCC_GAMMA_BODY, GammaCirc->data[j]);
    }
    fclose(fp);
  }

  /*initialize binary parameters structure*/
  XLAL_INIT_MEM(minBinaryTemplate);
  XLAL_INIT_MEM(maxBinaryTemplate);
  XLAL_INIT_MEM(thisBinaryTemplate);
  XLAL_INIT_MEM(binaryTemplateSpacings);
  /*fill in minbinaryOrbitParams*/
  XLALGPSSetREAL8( &minBinaryTemplate.tp, uvar.orbitTimeAsc);
  minBinaryTemplate.argp = 0.0;
  minBinaryTemplate.asini = uvar.orbitAsiniSec;
  minBinaryTemplate.ecc = 0.0;
  minBinaryTemplate.period = uvar.orbitPSec;
  minBinaryTemplate.fkdot[0] = uvar.fStart;
  /*fill in maxBinaryParams*/
  XLALGPSSetREAL8( &maxBinaryTemplate.tp, uvar.orbitTimeAsc + uvar.orbitTimeAscBand);
  maxBinaryTemplate.argp = 0.0;
  maxBinaryTemplate.asini = uvar.orbitAsiniSec + uvar.orbitAsiniSecBand;
  maxBinaryTemplate.ecc = 0.0;
  maxBinaryTemplate.period = uvar.orbitPSec;
  maxBinaryTemplate.fkdot[0] = uvar.fStart + uvar.fBand;
  /*fill in thisBinaryTemplate*/
  XLALGPSSetREAL8( &thisBinaryTemplate.tp, uvar.orbitTimeAsc + 0.5 * uvar.orbitTimeAscBand);
  thisBinaryTemplate.argp = 0.0;
  thisBinaryTemplate.asini = 0.5*(minBinaryTemplate.asini + maxBinaryTemplate.asini);
  thisBinaryTemplate.ecc = 0.0;
  thisBinaryTemplate.period =0.5*(minBinaryTemplate.period + maxBinaryTemplate.period);
  thisBinaryTemplate.fkdot[0]=0.5*(minBinaryTemplate.fkdot[0] + maxBinaryTemplate.fkdot[0]);

  /*Get metric diagonal components, also estimate sensitivity i.e. E[rho]/(h0)^2 (4.13)*/
  if ( (XLALCalculateLMXBCrossCorrDiagMetric(&estSens, &diagff, &diagaa, &diagTT, thisBinaryTemplate, GammaAve, sftPairs, sftIndices, inputSFTs, multiWeights /*, kappaValues*/)  != XLAL_SUCCESS ) ) {
    LogPrintf ( LOG_CRITICAL, "%s: XLALCalculateLMXBCrossCorrDiagMetric() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* spacing in frequency from diagff */ /* set spacings in new dopplerparams struct */
  if (XLALUserVarWasSet(&uvar.spacingF)) /* If spacing was given by CMD line, use it, else calculate spacing by mismatch*/
    binaryTemplateSpacings.fkdot[0] = uvar.spacingF;
  else
    binaryTemplateSpacings.fkdot[0] = sqrt(uvar.mismatchF / diagff);

  if (XLALUserVarWasSet(&uvar.spacingA))
    binaryTemplateSpacings.asini = uvar.spacingA;
  else
    binaryTemplateSpacings.asini = sqrt(uvar.mismatchA / diagaa);
  /* this is annoying: tp is a GPS time while we want a difference
     in time which should be just REAL8 */
  if (XLALUserVarWasSet(&uvar.spacingT))
    XLALGPSSetREAL8( &binaryTemplateSpacings.tp, uvar.spacingT);
  else
    XLALGPSSetREAL8( &binaryTemplateSpacings.tp, sqrt(uvar.mismatchT / diagTT));

  if (XLALUserVarWasSet(&uvar.spacingP))
    binaryTemplateSpacings.period = uvar.spacingP;
  else
    binaryTemplateSpacings.period = sqrt(uvar.mismatchP / diagpp);

  /* metric elements for eccentric case not considered? */

  UINT8 fCount = 0, aCount = 0, tCount = 0 , pCount = 0;
  const UINT8 fSpacingNum = floor( uvar.fBand / binaryTemplateSpacings.fkdot[0]);
  const UINT8 aSpacingNum = floor( uvar.orbitAsiniSecBand / binaryTemplateSpacings.asini);
  const UINT8 tSpacingNum = floor( uvar.orbitTimeAscBand / XLALGPSGetREAL8(&binaryTemplateSpacings.tp));
  const UINT8 pSpacingNum = floor( uvar.orbitPSecBand / binaryTemplateSpacings.period);

  /*reset minbinaryOrbitParams to shift the first point a factor so as to make the center of all seaching points centers at the center of searching band*/
  minBinaryTemplate.fkdot[0] = uvar.fStart + 0.5 * (uvar.fBand - fSpacingNum * binaryTemplateSpacings.fkdot[0]);
  minBinaryTemplate.asini = uvar.orbitAsiniSec + 0.5 * (uvar.orbitAsiniSecBand - aSpacingNum * binaryTemplateSpacings.asini);
  XLALGPSSetREAL8( &minBinaryTemplate.tp, uvar.orbitTimeAsc + 0.5 * (uvar.orbitTimeAscBand - tSpacingNum * XLALGPSGetREAL8(&binaryTemplateSpacings.tp)));
  minBinaryTemplate.period = uvar.orbitPSec + 0.5 * (uvar.orbitPSecBand - pSpacingNum * binaryTemplateSpacings.period);

  /* initialize the doppler scan struct which stores the current template information */
  XLALGPSSetREAL8(&dopplerpos.refTime, config.refTime);
  dopplerpos.Alpha = uvar.alphaRad;
  dopplerpos.Delta = uvar.deltaRad;
  dopplerpos.fkdot[0] = minBinaryTemplate.fkdot[0];
  /* set all spindowns to zero */
  for (k=1; k < PULSAR_MAX_SPINS; k++)
    dopplerpos.fkdot[k] = 0.0;
  dopplerpos.asini = minBinaryTemplate.asini;
  dopplerpos.period = minBinaryTemplate.period;
  dopplerpos.tp = minBinaryTemplate.tp;
  dopplerpos.ecc = minBinaryTemplate.ecc;
  dopplerpos.argp = minBinaryTemplate.argp;

  /* now set the initial values of binary parameters */
  /*  thisBinaryTemplate.asini = uvar.orbitAsiniSec;
  thisBinaryTemplate.period = uvar.orbitPSec;
  XLALGPSSetREAL8( &thisBinaryTemplate.tp, uvar.orbitTimeAsc);
  thisBinaryTemplate.ecc = 0.0;
  thisBinaryTemplate.argp = 0.0;*/
  /* copy to dopplerpos */

  /* Calculate SSB times (can do this once since search is currently only for one sky position, and binary doppler shift is added later) */
  MultiSSBtimes *multiSSBTimes = NULL;
  if ((multiSSBTimes = XLALGetMultiSSBtimes ( multiStates, skyPos, dopplerpos.refTime, SSBPREC_RELATIVISTICOPT )) == NULL){
    LogPrintf ( LOG_CRITICAL, "%s: XLALGetMultiSSBtimes() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* "New" general metric computation */
  /* For now hard-code circular parameter space */

  const DopplerCoordinateSystem coordSys = {
    .dim = 4,
    .coordIDs = { DOPPLERCOORD_FREQ,
		  DOPPLERCOORD_ASINI,
		  DOPPLERCOORD_TASC,
		  DOPPLERCOORD_PORB, },
  };

  REAL8VectorSequence *phaseDerivs = NULL;
  if ( ( XLALCalculateCrossCorrPhaseDerivatives ( &phaseDerivs, &thisBinaryTemplate, config.edat, sftIndices, multiSSBTimes, &coordSys )  != XLAL_SUCCESS ) ) {
    LogPrintf ( LOG_CRITICAL, "%s: XLALCalculateCrossCorrPhaseDerivatives() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* fill in metric and parameter offsets */
  gsl_matrix *g_ij = NULL;
  gsl_vector *eps_i = NULL;
  REAL8 sumGammaSq = 0;
  if ( ( XLALCalculateCrossCorrPhaseMetric ( &g_ij, &eps_i, &sumGammaSq, phaseDerivs, sftPairs, GammaAve, GammaCirc, &coordSys ) != XLAL_SUCCESS ) ) {
    LogPrintf ( LOG_CRITICAL, "%s: XLALCalculateCrossCorrPhaseMetric() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }
  XLALDestroyREAL8VectorSequence ( phaseDerivs );
  XLALDestroyREAL8Vector ( GammaCirc );

  if ((fp = fopen("gsldata.dat","w"))==NULL){
    LogPrintf ( LOG_CRITICAL, "Can't write in gsl matrix file");
    XLAL_ERROR( XLAL_EFUNC );
  }

  XLALfprintfGSLvector(fp, "%g", eps_i);
  XLALfprintfGSLmatrix(fp, "%g", g_ij);

  /* Allocate structure for binary doppler-shifting information */
  if ((multiBinaryTimes = XLALDuplicateMultiSSBtimes ( multiSSBTimes )) == NULL){
    LogPrintf ( LOG_CRITICAL, "%s: XLALDuplicateMultiSSBtimes() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  UINT8 numSFTs = sftIndices->length;
  if ((shiftedFreqs = XLALCreateREAL8Vector ( numSFTs ) ) == NULL){
    LogPrintf ( LOG_CRITICAL, "%s: XLALCreateREAL8Vector() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }
  if ((lowestBins = XLALCreateUINT4Vector ( numSFTs ) ) == NULL){
    LogPrintf ( LOG_CRITICAL, "%s: XLALCreateUINT4Vector() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  if ((expSignalPhases = XLALCreateCOMPLEX8Vector ( numSFTs ) ) == NULL){
    LogPrintf ( LOG_CRITICAL, "%s: XLALCreateREAL8Vector() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }
  if ((sincList = XLALCreateREAL8VectorSequence ( numSFTs, uvar.numBins ) ) == NULL){
    LogPrintf ( LOG_CRITICAL, "%s: XLALCreateREAL8VectorSequence() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* args should be : spacings, min and max doppler params */
  BOOLEAN firstPoint = TRUE; /* a boolean to help to search at the beginning point in parameter space, after the search it is set to be FALSE to end the loop*/
  if ( (XLALAddMultiBinaryTimes( &multiBinaryTimes, multiSSBTimes, &dopplerpos )  != XLAL_SUCCESS ) ) {
    LogPrintf ( LOG_CRITICAL, "%s: XLALAddMultiBinaryTimes() failed with errno=%d\n", __func__, xlalErrno );
    XLAL_ERROR( XLAL_EFUNC );
  } /*Need to apply additional doppler shifting before the loop, or the first point in parameter space will be lost and return a wrong SNR when fBand!=0*/

  while ( GetNextCrossCorrTemplate(&dopplerShiftFlag, &firstPoint, &dopplerpos, &binaryTemplateSpacings, &minBinaryTemplate, &maxBinaryTemplate, &fCount, &aCount, &tCount, &pCount, fSpacingNum, aSpacingNum, tSpacingNum, pSpacingNum) == 0)
    {
      /* do useful stuff here*/

      /* Apply additional Doppler shifting using current binary orbital parameters */
      /* Might want to be clever about checking whether we've changed the orbital parameters or only the frequency */
      if (dopplerShiftFlag == TRUE)
	{
	  if ( (XLALAddMultiBinaryTimes( &multiBinaryTimes, multiSSBTimes, &dopplerpos )  != XLAL_SUCCESS ) ) {
	    LogPrintf ( LOG_CRITICAL, "%s: XLALAddMultiBinaryTimes() failed with errno=%d\n", __func__, xlalErrno );
	    XLAL_ERROR( XLAL_EFUNC );
	  }
	}

      if ( (XLALGetDopplerShiftedFrequencyInfo( shiftedFreqs, lowestBins, expSignalPhases, sincList, uvar.numBins, &dopplerpos, sftIndices, inputSFTs, multiBinaryTimes, Tsft )  != XLAL_SUCCESS ) ) {
	LogPrintf ( LOG_CRITICAL, "%s: XLALGetDopplerShiftedFrequencyInfo() failed with errno=%d\n", __func__, xlalErrno );
	XLAL_ERROR( XLAL_EFUNC );
      }

      if ( (XLALCalculatePulsarCrossCorrStatistic( &ccStat, &evSquared, GammaAve, expSignalPhases, lowestBins, sincList, sftPairs, sftIndices, inputSFTs, multiWeights, uvar.numBins)  != XLAL_SUCCESS ) ) {
	LogPrintf ( LOG_CRITICAL, "%s: XLALCalculatePulsarCrossCorrStatistic() failed with errno=%d\n", __func__, xlalErrno );
	XLAL_ERROR( XLAL_EFUNC );
      }

      /* fill candidate struct and insert into toplist if necessary */
      thisCandidate.freq = dopplerpos.fkdot[0];
      thisCandidate.tp = XLALGPSGetREAL8( &dopplerpos.tp );
      thisCandidate.argp = dopplerpos.argp;
      thisCandidate.asini = dopplerpos.asini;
      thisCandidate.ecc = dopplerpos.ecc;
      thisCandidate.period = dopplerpos.period;
      thisCandidate.rho = ccStat;
      thisCandidate.evSquared = evSquared;
      thisCandidate.estSens = estSens;

      insert_into_crossCorrBinary_toplist(ccToplist, thisCandidate);

    } /* end while loop over templates */

  /* write candidates to file */
  sort_crossCorrBinary_toplist( ccToplist );
  /* add error checking */

  final_write_crossCorrBinary_toplist_to_file( ccToplist, uvar.toplistFilename, &checksum);

  REAL8 h0Sens = sqrt((10 / sqrt(estSens))); /*for a SNR=10 signal, the h0 we can detect*/

  XLALGPSTimeNow (&computingEndGPSTime); /*record the rough end time*/
  UINT4 computingTime = computingEndGPSTime.gpsSeconds - computingStartGPSTime.gpsSeconds;
  /* make a meta-data file*/
  if(XLALUserVarWasSet(&uvar.logFilename)){
    CHAR *CMDInputStr = XLALUserVarGetLog ( UVAR_LOGFMT_CFGFILE );
    if ((fp = fopen(uvar.logFilename,"w"))==NULL){
    LogPrintf ( LOG_CRITICAL, "Can't write in logfile");
    XLAL_ERROR( XLAL_EFUNC );
    }
    fprintf(fp, "[UserInput]\n\n");
    fprintf(fp, "%s\n", CMDInputStr);
    fprintf(fp, "[CalculatedValues]\n\n");
    fprintf(fp, "g_ff = %.9f\n", diagff );
    fprintf(fp, "g_aa = %.9f\n", diagaa );
    fprintf(fp, "g_TT = %.9f\n", diagTT );
    fprintf(fp, "FSpacing = %.9g\n", binaryTemplateSpacings.fkdot[0]);
    fprintf(fp, "ASpacing = %.9g\n", binaryTemplateSpacings.asini);
    fprintf(fp, "TSpacing = %.9g\n", XLALGPSGetREAL8(&binaryTemplateSpacings.tp));
    /* fprintf(fp, "PSpacing = %.9g\n", binaryTemplateSpacings.period );*/
    fprintf(fp, "TemplatenumF = %" LAL_UINT8_FORMAT "\n", (fSpacingNum + 1));
    fprintf(fp, "TemplatenumA = %" LAL_UINT8_FORMAT "\n", (aSpacingNum + 1));
    fprintf(fp, "TemplatenumT = %" LAL_UINT8_FORMAT "\n", (tSpacingNum + 1));
    fprintf(fp, "TemplatenumP = %" LAL_UINT8_FORMAT "\n", (pSpacingNum + 1));
    fprintf(fp, "TemplatenumTotal = %" LAL_UINT8_FORMAT "\n",(fSpacingNum + 1) * (aSpacingNum + 1) * (tSpacingNum + 1) * (pSpacingNum + 1));
    fprintf(fp, "Sens = %.9g\n", estSens);/*(E[rho]/h0^2)^2*/
    fprintf(fp, "h0_min_SNR10 = %.9g\n", h0Sens);/*for rho = 10 in our pipeline*/
    fprintf(fp, "startTime = %" LAL_INT4_FORMAT "\n", computingStartGPSTime.gpsSeconds );/*start time in GPS-time*/
    fprintf(fp, "endTime = %" LAL_INT4_FORMAT "\n", computingEndGPSTime.gpsSeconds );/*end time in GPS-time*/
    fprintf(fp, "computingTime = %" LAL_UINT4_FORMAT "\n", computingTime );/*total time in sec*/
    fprintf(fp, "SFTnum = %" LAL_UINT4_FORMAT "\n", sftIndices->length);/*total number of SFT*/
    fprintf(fp, "pairnum = %" LAL_UINT4_FORMAT "\n", sftPairs->length);/*total number of pair of SFT*/
    fprintf(fp, "Tsft = %.6g\n", Tsft);/*SFT duration*/
    fprintf(fp, "\n[Version]\n\n");
    fprintf(fp, "%s",  VCSInfoString);
    fclose(fp);
    XLALFree(CMDInputStr);
  }

  XLALFree(VCSInfoString);
  XLALDestroyCOMPLEX8Vector ( expSignalPhases );
  XLALDestroyUINT4Vector ( lowestBins );
  XLALDestroyREAL8Vector ( shiftedFreqs );
  XLALDestroyREAL8VectorSequence ( sincList );
  XLALDestroyMultiSSBtimes ( multiBinaryTimes );
  XLALDestroyMultiSSBtimes ( multiSSBTimes );
  XLALDestroyREAL8Vector ( GammaAve );
  XLALDestroySFTPairIndexList( sftPairs );
  XLALDestroySFTIndexList( sftIndices );
  XLALDestroyMultiAMCoeffs ( multiCoeffs );
  XLALDestroyMultiDetectorStateSeries ( multiStates );
  XLALDestroyMultiTimestamps ( multiTimes );
  XLALDestroyMultiNoiseWeights ( multiWeights );
  XLALDestroyMultiPSDVector ( multiPSDs );
  XLALDestroyMultiSFTVector ( inputSFTs );

  /* de-allocate memory for configuration variables */
  XLALDestroyConfigVars ( &config );

  /* de-allocate memory for user input variables */
  XLALDestroyUserVars();

  /* free toplist memory */
  free_crossCorr_toplist(&ccToplist);

  /* check memory leaks if we forgot to de-allocate anything */
  LALCheckMemoryLeaks();

  LogPrintf (LOG_CRITICAL, "End time\n");/*for debug convenience to record calculating time*/

  return 0;


} /* main */


/* initialize and register user variables */
int XLALInitUserVars (UserInput_t *uvar)
{

  /* initialize with some defaults */
  uvar->help = FALSE;
  uvar->maxLag = 0.0;
  uvar->inclAutoCorr = FALSE;
  uvar->fStart = 100.0;
  uvar->fBand = 0.1;
  /* uvar->fdotStart = 0.0; */
  /* uvar->fdotBand = 0.0; */
  uvar->alphaRad = 0.0;
  uvar->deltaRad = 0.0;
  uvar->refTime = 0.0;
  uvar->rngMedBlock = 50;
  uvar->numBins = 1;

  /* zero binary orbital parameters means not a binary */
  uvar->orbitAsiniSec = 0.0;
  uvar->orbitAsiniSecBand = 0.0;
  uvar->orbitPSec = 0.0;
  uvar->orbitPSecBand = 0.0;
  uvar->orbitTimeAsc = 0;
  uvar->orbitTimeAscBand = 0;

  /*default mismatch values */
  /* set to 0.1 by default -- for no real reason */
  /* make 0.1 a macro? */
  uvar->mismatchF = 0.1;
  uvar->mismatchA = 0.1;
  uvar->mismatchT = 0.1;
  uvar->mismatchP = 0.1;

  uvar->ephemEarth = XLALStringDuplicate("earth00-19-DE405.dat.gz");
  uvar->ephemSun = XLALStringDuplicate("sun00-19-DE405.dat.gz");

  uvar->sftLocation = XLALCalloc(1, MAXFILENAMELENGTH+1);

  /* initialize number of candidates in toplist -- default is just to return the single best candidate */
  uvar->numCand = 1;
  uvar->toplistFilename = XLALStringDuplicate("toplist_crosscorr.dat");
  uvar->version = FALSE;

  /* register  user-variables */
  XLALregBOOLUserStruct  ( help, 	   'h',  UVAR_HELP, "Print this message");
  XLALregINTUserStruct   ( startTime,       0,  UVAR_REQUIRED, "Desired start time of analysis in GPS seconds");
  XLALregINTUserStruct   ( endTime,         0,  UVAR_REQUIRED, "Desired end time of analysis in GPS seconds");
  XLALregREALUserStruct  ( maxLag,          0,  UVAR_OPTIONAL, "Maximum lag time in seconds between SFTs in correlation");
  XLALregBOOLUserStruct  ( inclAutoCorr,    0,  UVAR_OPTIONAL, "Include auto-correlation terms (an SFT with itself)");
  XLALregREALUserStruct  ( fStart,          0,  UVAR_OPTIONAL, "Start frequency in Hz");
  XLALregREALUserStruct  ( fBand,           0,  UVAR_OPTIONAL, "Frequency band to search over in Hz ");
  /* XLALregREALUserStruct  ( fdotStart,     0,  UVAR_OPTIONAL, "Start value of spindown in Hz/s"); */
  /* XLALregREALUserStruct  ( fdotBand,      0,  UVAR_OPTIONAL, "Band for spindown values in Hz/s"); */
  XLALregREALUserStruct  ( alphaRad,        0,  UVAR_OPTIONAL, "Right ascension for directed search (radians)");
  XLALregREALUserStruct  ( deltaRad,        0,  UVAR_OPTIONAL, "Declination for directed search (radians)");
  XLALregREALUserStruct  ( refTime,         0,  UVAR_OPTIONAL, "SSB reference time for pulsar-parameters [Default: midPoint]");
  XLALregREALUserStruct  ( orbitAsiniSec,   0,  UVAR_OPTIONAL, "Start of search band for projected semimajor axis (seconds) [0 means not a binary]");
  XLALregREALUserStruct  ( orbitAsiniSecBand, 0,  UVAR_OPTIONAL, "Width of search band for projected semimajor axis (seconds)");
  XLALregREALUserStruct  ( orbitPSec,       0,  UVAR_OPTIONAL, "Binary orbital period (seconds) [0 means not a binary]");
  XLALregREALUserStruct  ( orbitPSecBand,       0,  UVAR_OPTIONAL, "Band for binary orbital period (seconds) ");
  XLALregREALUserStruct  ( orbitTimeAsc,    0,  UVAR_OPTIONAL, "Start of orbital time-of-ascension band in GPS seconds");
  XLALregREALUserStruct  ( orbitTimeAscBand, 0,  UVAR_OPTIONAL, "Width of orbital time-of-ascension band (seconds)");
  XLALregSTRINGUserStruct( ephemEarth,      0,  UVAR_OPTIONAL, "Earth ephemeris file to use");
  XLALregSTRINGUserStruct( ephemSun,        0,  UVAR_OPTIONAL, "Sun ephemeris file to use");
  XLALregSTRINGUserStruct( sftLocation,     0,  UVAR_REQUIRED, "Filename pattern for locating SFT data");
  XLALregINTUserStruct   ( rngMedBlock,     0,  UVAR_OPTIONAL, "Running median block size for PSD estimation");
  XLALregINTUserStruct   ( numBins,         0,  UVAR_OPTIONAL, "Number of frequency bins to include in calculation");
  XLALregREALUserStruct  ( mismatchF,       0,  UVAR_OPTIONAL, "Desired mismatch for frequency spacing");
  XLALregREALUserStruct  ( mismatchA,       0,  UVAR_OPTIONAL, "Desired mismatch for asini spacing");
  XLALregREALUserStruct  ( mismatchT,       0,  UVAR_OPTIONAL, "Desired mismatch for periapse passage time spacing");
  XLALregREALUserStruct  ( mismatchP,       0,  UVAR_OPTIONAL, "Desired mismatch for period spacing");
  XLALregREALUserStruct  ( spacingF,       0,  UVAR_OPTIONAL, "Desired frequency spacing");
  XLALregREALUserStruct  ( spacingA,       0,  UVAR_OPTIONAL, "Desired asini spacing");
  XLALregREALUserStruct  ( spacingT,       0,  UVAR_OPTIONAL, "Desired periapse passage time spacing");
  XLALregREALUserStruct  ( spacingP,       0,  UVAR_OPTIONAL, "Desired period spacing");
  XLALregINTUserStruct   ( numCand,         0,  UVAR_OPTIONAL, "Number of candidates to keep in toplist");
  XLALregSTRINGUserStruct( pairListInputFilename, 0,  UVAR_OPTIONAL, "Name of file from which to read list of SFT pairs");
  XLALregSTRINGUserStruct( pairListOutputFilename, 0,  UVAR_OPTIONAL, "Name of file to which to write list of SFT pairs");
  XLALregSTRINGUserStruct( sftListOutputFilename, 0,  UVAR_OPTIONAL, "Name of file to which to write list of SFTs (for sanity checks)");
  XLALregSTRINGUserStruct( sftListInputFilename, 0,  UVAR_OPTIONAL, "Name of file to which to read in list of SFTs (for sanity checks)");
  XLALregSTRINGUserStruct( gammaAveOutputFilename, 0,  UVAR_OPTIONAL, "Name of file to which to write aa+bb weights (for e.g., false alarm estimation)");
  XLALregSTRINGUserStruct( gammaCircOutputFilename, 0,  UVAR_OPTIONAL, "Name of file to which to write ab-ba weights (for e.g., systematic error)");
  XLALregSTRINGUserStruct( toplistFilename, 0,  UVAR_OPTIONAL, "Output filename containing candidates in toplist");
  XLALregSTRINGUserStruct( logFilename, 0,  UVAR_OPTIONAL, "Output a meta-data file for the search");
  XLALregBOOLUserStruct  ( version, 	   'V',  UVAR_SPECIAL, "Output version(VCS) information");
  if ( xlalErrno ) {
    XLALPrintError ("%s: user variable initialization failed with errno = %d.\n", __func__, xlalErrno );
    XLAL_ERROR ( XLAL_EFUNC );
  }

  return XLAL_SUCCESS;
}
REAL8VectorSequence * XLALASCIIFileReadColumns( INT4 ncol, const char *fname )
{
  char line[LINE_MAX];
  REAL8VectorSequence *data;
  int nline;
  int nrow;
  int row;
  int col;
  FILE *fp;

  if ( ncol < 0 )
    XLAL_ERROR_NULL( XLAL_EINVAL );

  /* count rows */
  nrow = XLALASCIIFileCountRows( fname );
  if ( nrow < 0 )
    XLAL_ERROR_NULL( XLAL_EFUNC );

  /* allocate memory for data */
  /* column 0 will contain line number of input file */
  data = XLALCreateREAL8VectorSequence( nrow, ncol + 1 );
  if ( ! data )
    XLAL_ERROR_NULL( XLAL_EFUNC );

  /* open file */
  fp = fopen( fname, "r" );
  if ( ! fp )
  {
    XLALDestroyREAL8VectorSequence( data );
    XLAL_ERROR_NULL( XLAL_EIO );
  }

  nline = 0;
  row = 0;
  while ( row < nrow )
    if ( fgets( line, sizeof( line ), fp ) )
    {
      char *p = line;
      ++nline;
      if ( strlen( line ) >= sizeof( line ) - 1 )
      {
        XLALPrintError( "XLAL Error - %s: line %d too long\n\tfile: %s\n", __func__, nline, fname );
        XLALDestroyREAL8VectorSequence( data );
        fclose( fp );
        XLAL_ERROR_NULL( XLAL_EBADLEN );
      }
      if ( line[0] == '#' || line[0] == '%' )
        continue;
      data->data[row*data->vectorLength] = nline;
      for ( col = 1; col <= ncol; ++col )
      {
        char *endp;
        REAL8 val;
        val = strtod( p, &endp );
        if ( p == endp ) /* no conversion */
        {
          XLALPrintError( "XLAL Error - %s: unable to parse line %d\n\tfile: %s\n", __func__, nline, fname );
          XLALDestroyREAL8VectorSequence( data );
          fclose( fp );
          XLAL_ERROR_NULL( XLAL_EFAILED );
        }
        if ( isnan( val ) )
        {
          XLALPrintWarning( "XLAL Warning - %s: invalid data (nan) in column %d of line %d\n\tfile: %s\n", __func__, col, nline, fname );
          val = 0;
        }
        if ( isinf( val ) )
        {
          XLALPrintWarning( "XLAL Warning - %s: invalid data (inf) in column %d of line %d\n\tfile: %s\n", __func__, col, nline, fname );
          val = 0;
        }
        data->data[row*data->vectorLength + col] = val;
        p = endp;
      }
      ++row;
    }
    else
    {
      XLALDestroyREAL8VectorSequence( data );
      XLAL_ERROR_NULL( XLAL_EIO );
    }

  fclose( fp );
  return data;
}
/**
 * The main workhorse function for performing the ringdown attachment for EOB
 * models EOBNRv2 and SEOBNRv1. This is the function which gets called by the
 * code generating the full IMR waveform once generation of the inspiral part
 * has been completed.
 * The ringdown is attached using the hybrid comb matching detailed in
 * The method is describe in Sec. II C of Pan et al. PRD 84, 124052 (2011),
 * specifically Eqs. 30 - 32.. Further details of the
 * implementation of the found in the DCC document T1100433.
 * In SEOBNRv1, the last physical overtone is replace by a pseudoQNM. See
 * Taracchini et al. PRD 86, 024011 (2012) for details.
 * STEP 1) Get mass and spin of the final black hole and the complex ringdown frequencies
 * STEP 2) Based on least-damped-mode decay time, allocate memory for rigndown waveform
 * STEP 3) Get values and derivatives of inspiral waveforms at matching comb points
 * STEP 4) Solve QNM coefficients and generate ringdown waveforms
 * STEP 5) Stitch inspiral and ringdown waveoforms
 */
static INT4 XLALSimIMREOBHybridAttachRingdown(
  REAL8Vector *signal1,    /**<< OUTPUT, Real of inspiral waveform to which we attach ringdown */
  REAL8Vector *signal2,    /**<< OUTPUT, Imag of inspiral waveform to which we attach ringdown */
  const INT4   l,          /**<< Current mode l */
  const INT4   m,          /**<< Current mode m */
  const REAL8  dt,         /**<< Sample time step (in seconds) */
  const REAL8  mass1,      /**<< First component mass (in Solar masses) */
  const REAL8  mass2,      /**<< Second component mass (in Solar masses) */
  const REAL8  spin1x,     /**<<The spin of the first object; only needed for spin waveforms */
  const REAL8  spin1y,     /**<<The spin of the first object; only needed for spin waveforms */
  const REAL8  spin1z,     /**<<The spin of the first object; only needed for spin waveforms */
  const REAL8  spin2x,     /**<<The spin of the second object; only needed for spin waveforms */
  const REAL8  spin2y,     /**<<The spin of the second object; only needed for spin waveforms */
  const REAL8  spin2z,     /**<<The spin of the second object; only needed for spin waveforms */
  REAL8Vector *timeVec,    /**<< Vector containing the time values */
  REAL8Vector *matchrange, /**<< Time values chosen as points for performing comb matching */
  Approximant  approximant /**<<The waveform approximant being used */
  )
{

      COMPLEX16Vector *modefreqs;
      UINT4 Nrdwave;
      UINT4 j;

      UINT4 nmodes;
      REAL8Vector		*rdwave1;
      REAL8Vector		*rdwave2;
      REAL8Vector		*rinspwave;
      REAL8Vector		*dinspwave;
      REAL8Vector		*ddinspwave;
      REAL8VectorSequence	*inspwaves1;
      REAL8VectorSequence	*inspwaves2;
      REAL8 eta, a, NRPeakOmega22; /* To generate pQNM frequency */
      REAL8 mTot; /* In geometric units */
      REAL8 spin1[3] = { spin1x, spin1y, spin1z };
      REAL8 spin2[3] = { spin2x, spin2y, spin2z };
      REAL8 finalMass, finalSpin;

      mTot  = (mass1 + mass2) * LAL_MTSUN_SI;
      eta       = mass1 * mass2 / ( (mass1 + mass2) * (mass1 + mass2) );

      /*
       * STEP 1) Get mass and spin of the final black hole and the complex ringdown frequencies
       */

      /* Create memory for the QNM frequencies */
      nmodes = 8;
      modefreqs = XLALCreateCOMPLEX16Vector( nmodes );
      if ( !modefreqs )
      {
        XLAL_ERROR( XLAL_ENOMEM );
      }

      if ( XLALSimIMREOBGenerateQNMFreqV2( modefreqs, mass1, mass2, spin1, spin2, l, m, nmodes, approximant ) == XLAL_FAILURE )
      {
        XLALDestroyCOMPLEX16Vector( modefreqs );
        XLAL_ERROR( XLAL_EFUNC );
      }

      /* Call XLALSimIMREOBFinalMassSpin() to get mass and spin of the final black hole */
      if ( XLALSimIMREOBFinalMassSpin(&finalMass, &finalSpin, mass1, mass2, spin1, spin2, approximant) == XLAL_FAILURE )
      {
        XLAL_ERROR( XLAL_EFUNC );
      }

      if ( approximant == SEOBNRv1 )
      {
          /* Replace the last QNM with pQNM */
          /* We assume aligned/antialigned spins here */
          a  = (spin1[2] + spin2[2]) / 2. * (1.0 - 2.0 * eta) + (spin1[2] - spin2[2]) / 2. * (mass1 - mass2) / (mass1 + mass2);
          NRPeakOmega22 = GetNRSpinPeakOmega( l, m, eta, a ) / mTot;
          /*printf("a and NRomega in QNM freq: %.16e %.16e %.16e %.16e %.16e\n",spin1[2],spin2[2],
                 mTot/LAL_MTSUN_SI,a,NRPeakOmega22*mTot);*/
          modefreqs->data[7] = (NRPeakOmega22/finalMass + creal(modefreqs->data[0])) / 2.;
          modefreqs->data[7] += I * 10./3. * cimag(modefreqs->data[0]);
      }

      /*for (j = 0; j < nmodes; j++)
      {
        printf("QNM frequencies: %d %d %d %e %e\n",l,m,j,modefreqs->data[j].re*mTot,1./modefreqs->data[j].im/mTot);
      }*/

      /* Ringdown signal length: 10 times the decay time of the n=0 mode */
      Nrdwave = (INT4) (EOB_RD_EFOLDS / cimag(modefreqs->data[0]) / dt);

      /* Check the value of attpos, to prevent memory access problems later */
      if ( matchrange->data[0] * mTot / dt < 5 || matchrange->data[1]*mTot/dt > matchrange->data[2] *mTot/dt - 2 )
      {
        XLALPrintError( "More inspiral points needed for ringdown matching.\n" );
        //printf("%.16e,%.16e,%.16e\n",matchrange->data[0] * mTot / dt, matchrange->data[1]*mTot/dt, matchrange->data[2] *mTot/dt - 2);
        XLALDestroyCOMPLEX16Vector( modefreqs );
        XLAL_ERROR( XLAL_EFAILED );
      }

      /*
       * STEP 2) Based on least-damped-mode decay time, allocate memory for rigndown waveform
       */

      /* Create memory for the ring-down and full waveforms, and derivatives of inspirals */

      rdwave1 = XLALCreateREAL8Vector( Nrdwave );
      rdwave2 = XLALCreateREAL8Vector( Nrdwave );
      rinspwave = XLALCreateREAL8Vector( 6 );
      dinspwave = XLALCreateREAL8Vector( 6 );
      ddinspwave = XLALCreateREAL8Vector( 6 );
      inspwaves1 = XLALCreateREAL8VectorSequence( 3, 6 );
      inspwaves2 = XLALCreateREAL8VectorSequence( 3, 6 );

      /* Check memory was allocated */
      if ( !rdwave1 || !rdwave2 || !rinspwave || !dinspwave 
	   || !ddinspwave || !inspwaves1 || !inspwaves2 )
      {
        XLALDestroyCOMPLEX16Vector( modefreqs );
        if (rdwave1)    XLALDestroyREAL8Vector( rdwave1 );
        if (rdwave2)    XLALDestroyREAL8Vector( rdwave2 );
        if (rinspwave)  XLALDestroyREAL8Vector( rinspwave );
        if (dinspwave)  XLALDestroyREAL8Vector( dinspwave );
        if (ddinspwave) XLALDestroyREAL8Vector( ddinspwave );
        if (inspwaves1) XLALDestroyREAL8VectorSequence( inspwaves1 );
        if (inspwaves2) XLALDestroyREAL8VectorSequence( inspwaves2 );
        XLAL_ERROR( XLAL_ENOMEM );
      }

      memset( rdwave1->data, 0, rdwave1->length * sizeof( REAL8 ) );
      memset( rdwave2->data, 0, rdwave2->length * sizeof( REAL8 ) );

      /*
       * STEP 3) Get values and derivatives of inspiral waveforms at matching comb points
       */

      /* Generate derivatives of the last part of inspiral waves */
      /* Get derivatives of signal1 */
      if ( XLALGenerateHybridWaveDerivatives( rinspwave, dinspwave, ddinspwave, timeVec, signal1, 
			matchrange, dt, mass1, mass2 ) == XLAL_FAILURE )
      {
        XLALDestroyCOMPLEX16Vector( modefreqs );
        XLALDestroyREAL8Vector( rdwave1 );
        XLALDestroyREAL8Vector( rdwave2 );
        XLALDestroyREAL8Vector( rinspwave );
        XLALDestroyREAL8Vector( dinspwave );
        XLALDestroyREAL8Vector( ddinspwave );
        XLALDestroyREAL8VectorSequence( inspwaves1 );
        XLALDestroyREAL8VectorSequence( inspwaves2 );
        XLAL_ERROR( XLAL_EFUNC );
      }
      for (j = 0; j < 6; j++)
      {
	    inspwaves1->data[j] = rinspwave->data[j];
	    inspwaves1->data[j + 6] = dinspwave->data[j];
	    inspwaves1->data[j + 12] = ddinspwave->data[j];
      }

      /* Get derivatives of signal2 */
      if ( XLALGenerateHybridWaveDerivatives( rinspwave, dinspwave, ddinspwave, timeVec, signal2, 
			matchrange, dt, mass1, mass2 ) == XLAL_FAILURE )
      {
        XLALDestroyCOMPLEX16Vector( modefreqs );
        XLALDestroyREAL8Vector( rdwave1 );
        XLALDestroyREAL8Vector( rdwave2 );
        XLALDestroyREAL8Vector( rinspwave );
        XLALDestroyREAL8Vector( dinspwave );
        XLALDestroyREAL8Vector( ddinspwave );
        XLALDestroyREAL8VectorSequence( inspwaves1 );
        XLALDestroyREAL8VectorSequence( inspwaves2 );
        XLAL_ERROR( XLAL_EFUNC );
      }
      for (j = 0; j < 6; j++)
      {
	    inspwaves2->data[j] = rinspwave->data[j];
	    inspwaves2->data[j + 6] = dinspwave->data[j];
	    inspwaves2->data[j + 12] = ddinspwave->data[j];
      }


      /*
       * STEP 4) Solve QNM coefficients and generate ringdown waveforms
       */

      /* Generate ring-down waveforms */
      if ( XLALSimIMREOBHybridRingdownWave( rdwave1, rdwave2, dt, mass1, mass2, inspwaves1, inspwaves2,
			  modefreqs, matchrange ) == XLAL_FAILURE )
      {
        XLALDestroyCOMPLEX16Vector( modefreqs );
        XLALDestroyREAL8Vector( rdwave1 );
        XLALDestroyREAL8Vector( rdwave2 );
        XLALDestroyREAL8Vector( rinspwave );
        XLALDestroyREAL8Vector( dinspwave );
        XLALDestroyREAL8Vector( ddinspwave );
        XLALDestroyREAL8VectorSequence( inspwaves1 );
        XLALDestroyREAL8VectorSequence( inspwaves2 );
        XLAL_ERROR( XLAL_EFUNC );
      }

      /*
       * STEP 5) Stitch inspiral and ringdown waveoforms
       */

      /* Generate full waveforms, by stitching inspiral and ring-down waveforms */
      UINT4 attachIdx = matchrange->data[1] * mTot / dt;
      for (j = 1; j < Nrdwave; ++j)
      {
	    signal1->data[j + attachIdx] = rdwave1->data[j];
	    signal2->data[j + attachIdx] = rdwave2->data[j];
      }

      memset( signal1->data+Nrdwave+attachIdx, 0, (signal1->length - Nrdwave - attachIdx)*sizeof(REAL8) );
      memset( signal2->data+Nrdwave+attachIdx, 0, (signal2->length - Nrdwave - attachIdx)*sizeof(REAL8) );

      /* Free memory */
      XLALDestroyCOMPLEX16Vector( modefreqs );
      XLALDestroyREAL8Vector( rdwave1 );
      XLALDestroyREAL8Vector( rdwave2 );
      XLALDestroyREAL8Vector( rinspwave );
      XLALDestroyREAL8Vector( dinspwave );
      XLALDestroyREAL8Vector( ddinspwave );
      XLALDestroyREAL8VectorSequence( inspwaves1 );
      XLALDestroyREAL8VectorSequence( inspwaves2 );

      return XLAL_SUCCESS;
}
/**
 * Helper funxtion to copy velocity, time and position vectors out of the
 * multi-detector state series
 */
void LALGetMultiDetectorVelTimePos(LALStatus                *status,
				   REAL8VectorSequence      **outVel,
				   REAL8VectorSequence      **outPos,
				   LIGOTimeGPSVector        **outTime,
				   MultiDetectorStateSeries *in)
{

  UINT4 numifo, len, numsft, iIFO, iSFT, j;
  REAL8VectorSequence *velV = NULL;
  REAL8VectorSequence *posV = NULL;
  LIGOTimeGPSVector *timeV = NULL;

  INITSTATUS(status);
  ATTATCHSTATUSPTR (status);

  ASSERT (in, status, DETECTORSTATES_ENULL, DETECTORSTATES_MSGENULL);
  ASSERT (in->length > 0, status, DETECTORSTATES_ENULL, DETECTORSTATES_MSGENULL);

  ASSERT (*outVel == NULL, status, DETECTORSTATES_ENONULL, DETECTORSTATES_MSGENONULL);
  ASSERT (*outPos == NULL, status, DETECTORSTATES_ENONULL, DETECTORSTATES_MSGENONULL);
  ASSERT (*outTime == NULL, status, DETECTORSTATES_ENONULL, DETECTORSTATES_MSGENONULL);

  /* number of ifos */
  numifo = in->length;

  len = 0;
  /* calculate total number of data points */
  for (j = 0, iIFO = 0; iIFO < numifo; iIFO++ ) {
    ASSERT (in->data[iIFO], status, DETECTORSTATES_ENULL, DETECTORSTATES_MSGENULL);
    len += in->data[iIFO]->length;
  }

  /* allocate memory for vectors */
  velV =  XLALCreateREAL8VectorSequence ( len, 3 );
  posV =  XLALCreateREAL8VectorSequence ( len, 3 );
  TRY (LALCreateTimestampVector ( status->statusPtr, &timeV,  len), status);

  /* copy the timestamps, weights, and velocity vector */
  for (j = 0, iIFO = 0; iIFO < numifo; iIFO++ ) {

    ASSERT (in->data[iIFO], status, DETECTORSTATES_ENULL, DETECTORSTATES_MSGENULL);

    numsft = in->data[iIFO]->length;
    for ( iSFT = 0; iSFT < numsft; iSFT++, j++) {

      velV->data[3*j] = in->data[iIFO]->data[iSFT].vDetector[0];
      velV->data[3*j+1] = in->data[iIFO]->data[iSFT].vDetector[1];
      velV->data[3*j+2] = in->data[iIFO]->data[iSFT].vDetector[2];

      posV->data[3*j] = in->data[iIFO]->data[iSFT].rDetector[0];
      posV->data[3*j+1] = in->data[iIFO]->data[iSFT].rDetector[1];
      posV->data[3*j+2] = in->data[iIFO]->data[iSFT].rDetector[2];

      /* mid time of sfts */
      timeV->data[j] = in->data[iIFO]->data[iSFT].tGPS;

    } /* loop over SFTs */

  } /* loop over IFOs */

  *outVel = velV;
  *outPos = posV;
  *outTime = timeV;

  DETATCHSTATUSPTR (status);

  /* normal exit */
  RETURN (status);
}
Beispiel #5
0
void LALReadNRWave_raw_real8(LALStatus *status,	/**< pointer to LALStatus structure */
		       REAL8TimeVectorSeries **out, /**< [out] output time series for h+ and hx */
		       const CHAR  *filename        /**< [in] File containing numrel waveform */)
{

  UINT4 length, k, r;
  REAL8TimeVectorSeries *ret=NULL;
  REAL8VectorSequence *data=NULL;
  REAL8Vector *timeVec=NULL;
  LALParsedDataFile *cfgdata=NULL;
  REAL8 tmp1, tmp2, tmp3;

  INITSTATUS(status);
  ATTATCHSTATUSPTR (status);

  /* some consistency checks */
  ASSERT (filename != NULL, status, NRWAVEIO_ENULL, NRWAVEIO_MSGENULL );
  ASSERT ( out != NULL, status, NRWAVEIO_ENULL, NRWAVEIO_MSGENULL );
  ASSERT ( *out == NULL, status, NRWAVEIO_ENONULL, NRWAVEIO_MSGENONULL );

  if ( XLALParseDataFile ( &cfgdata, filename ) != XLAL_SUCCESS ) {
    ABORT( status, NRWAVEIO_EFILE, NRWAVEIO_MSGEFILE );
  }
  length = cfgdata->lines->nTokens; /*number of data points */


  /* allocate memory */
  ret = LALCalloc(1, sizeof(*ret));
  if (!ret) {
    ABORT( status, NRWAVEIO_ENOMEM, NRWAVEIO_MSGENOMEM );
  }
  strcpy(ret->name,filename);
  ret->f0 = 0;

  data =  XLALCreateREAL8VectorSequence (2, length);
  if (!data) {
    ABORT( status, NRWAVEIO_ENOMEM, NRWAVEIO_MSGENOMEM );
  }

  timeVec = XLALCreateREAL8Vector (length);
  if (!timeVec) {
    ABORT( status, NRWAVEIO_ENOMEM, NRWAVEIO_MSGENOMEM );
  }

  /* now get the data */
  for (k = 0; k < length; k++) {
    r = sscanf(cfgdata->lines->tokens[k], "%lf%lf%lf", &tmp1, &tmp2, &tmp3);

    /* Check the data file format */
    if ( r != 3) {
      /* there must be exactly 3 data entries -- time, h+, hx */
      ABORT( status, NRWAVEIO_EFORMAT, NRWAVEIO_MSGEFORMAT );
    }

    timeVec->data[k] = tmp1;
    data->data[k] = tmp2;
    data->data[data->vectorLength + k] = tmp3;

  }

  /*  scale time */
  ret->deltaT = timeVec->data[1] - timeVec->data[0];

  /* might also want to go through timeVec to make sure it is evenly spaced */


  ret->data = data;
  (*out) = ret;

  XLALDestroyREAL8Vector (timeVec);
  XLALDestroyParsedDataFile (cfgdata);

  DETATCHSTATUSPTR(status);
  RETURN(status);

} /* LALReadNRWave() */