/** * Compute transient-CW Bayes-factor B_SG = P(x|HypS)/P(x|HypG) (where HypG = Gaussian noise hypothesis), * marginalized over start-time and timescale of transient CW signal, using given type and parameters * of transient window range. * * Note: this function is a C-implemention, partially based-on/inspired-by Stefanos Giampanis' * original matlab implementation of this search function. * * Note2: if window->type == none, uses a single rectangular window covering all the data. */ REAL8 XLALComputeTransientBstat ( transientWindowRange_t windowRange, /**< [in] type and parameters specifying transient window range */ const transientFstatMap_t *FstatMap /**< [in] pre-computed transient-Fstat map F_mn over {t0, tau} ranges */ ) { /* ----- check input consistency */ if ( !FstatMap || !FstatMap->F_mn ) { XLALPrintError ("%s: invalid NULL input 'FstatMap' or 'FstatMap->F_mn'\n", __func__ ); XLAL_ERROR ( XLAL_EINVAL ); } if ( windowRange.type >= TRANSIENT_LAST ) { XLALPrintError ("%s: unknown window-type (%d) passes as input. Allowed are [0,%d].\n", __func__, windowRange.type, TRANSIENT_LAST-1); XLAL_ERROR ( XLAL_EINVAL ); } /* ----- step through F_mn array subtract maxF and sum e^{F_mn - maxF}*/ /* * The maximum-likelihood Fmax is globally subtracted from F_mn, and stored separatedly in the struct, because in most * expressions it is numerically more robust to compute e^(F_mn - Fmax), which at worst can underflow, while * e^F_mn can overflow (for F>~700). The constant offset e^Fmax is irrelevant for posteriors (normalization constant), or * can be handled separately, eg by computing log(B) = Fmax + log(sum e^(Fmn-Fmax)) for the Bayes-factor. */ UINT4 N_t0Range = FstatMap->F_mn->size1; UINT4 N_tauRange = FstatMap->F_mn->size2; UINT4 m, n; REAL8 sum_eB = 0; for ( m=0; m < N_t0Range; m ++ ) { for ( n=0; n < N_tauRange; n ++ ) { REAL8 DeltaF = FstatMap->maxF - gsl_matrix_get ( FstatMap->F_mn, m, n ); // always >= 0, exactly ==0 at {m,n}_max //sum_eB += exp ( - DeltaF ); sum_eB += XLALFastNegExp ( DeltaF ); } /* for n < N_tauRange */ } /* for m < N_t0Range */ /* combine this to final log(Bstat) result with proper normalization (assuming rhohMax=1) : */ REAL8 logBhat = FstatMap->maxF + log ( sum_eB ); // unnormalized Bhat REAL8 normBh = 70.0 / ( N_t0Range * N_tauRange ); // normalization factor assuming rhohMax=1 /* final normalized Bayes factor, assuming rhohMax=1 */ /* NOTE: correct this for different rhohMax by adding "- 4 * log(rhohMax)" to logB*/ REAL8 logBstat = log ( normBh ) + logBhat; /* - 4.0 * log ( rhohMax ) */ // printf ( "\n\nlogBhat = %g, normBh = %g, log(normBh) = %g\nN_t0Range = %d, N_tauRange=%d\n\n", logBhat, normBh, log(normBh), N_t0Range, N_tauRange ); /* free mem */ XLALDestroyExpLUT(); /* ----- return ----- */ return logBstat; } /* XLALComputeTransientBstat() */
/** * MAIN function * Generates samples of B-stat and F-stat according to their pdfs for given signal-params. */ int main(int argc,char *argv[]) { UserInput_t XLAL_INIT_DECL(uvar); ConfigVariables XLAL_INIT_DECL(cfg); vrbflg = 1; /* verbose error-messages */ /* turn off default GSL error handler */ gsl_set_error_handler_off (); /* ----- register and read all user-variables ----- */ if ( XLALInitUserVars( &uvar ) != XLAL_SUCCESS ) { LogPrintf ( LOG_CRITICAL, "%s: XLALInitUserVars() failed with errno=%d\n", __func__, xlalErrno ); return 1; } /* do ALL cmdline and cfgfile handling */ BOOLEAN should_exit = 0; if ( XLALUserVarReadAllInput ( &should_exit, argc, argv ) != XLAL_SUCCESS ) { LogPrintf ( LOG_CRITICAL, "%s: XLALUserVarReadAllInput() failed with errno=%d\n", __func__, xlalErrno ); return 1; } if ( should_exit ) return EXIT_FAILURE; if ( uvar.version ) { /* output verbose VCS version string if requested */ CHAR *vcs; if ( (vcs = XLALGetVersionString (lalDebugLevel)) == NULL ) { LogPrintf ( LOG_CRITICAL, "%s:XLALGetVersionString(%d) failed with errno=%d.\n", __func__, lalDebugLevel, xlalErrno ); return 1; } printf ( "%s\n", vcs ); XLALFree ( vcs ); return 0; } /* ---------- Initialize code-setup ---------- */ if ( XLALInitCode( &cfg, &uvar ) != XLAL_SUCCESS ) { LogPrintf (LOG_CRITICAL, "%s: XLALInitCode() failed with error = %d\n", __func__, xlalErrno ); XLAL_ERROR ( XLAL_EFUNC ); } /* compare IFO name for line injection with IFO list, find the corresponding index, or throw an error if not found */ UINT4 numDetectors = cfg.multiDetStates->length; INT4 lineX = -1; if ( uvar.lineIFO ) { for ( UINT4 X=0; X < numDetectors; X++ ) { if ( strcmp( uvar.lineIFO, uvar.IFOs->data[X] ) == 0 ) lineX = X; } if ( lineX == -1 ) { XLALPrintError ("\nError in function %s, line %d : Could not match detector ID \"%s\" for line injection to any detector.\n\n", __func__, __LINE__, uvar.lineIFO); XLAL_ERROR ( XLAL_EFAILED ); } } /* ----- prepare stats output ----- */ FILE *fpStats = NULL; if ( uvar.outputStats ) { if ( (fpStats = fopen (uvar.outputStats, "wb")) == NULL) { LogPrintf (LOG_CRITICAL, "Error opening file '%s' for writing..\n\n", uvar.outputStats ); XLAL_ERROR ( XLAL_EIO ); } fprintf (fpStats, "%s", cfg.logString ); /* write search log comment */ if ( write_BSGL_candidate_to_fp ( fpStats, NULL, uvar.IFOs, NULL, uvar.computeBSGL ) != XLAL_SUCCESS ) { /* write header-line comment */ XLAL_ERROR ( XLAL_EFUNC ); } } /* if outputStats */ /* ----- prepare injection params output ----- */ FILE *fpInjParams = NULL; if ( uvar.outputInjParams ) { if ( (fpInjParams = fopen (uvar.outputInjParams, "wb")) == NULL) { LogPrintf (LOG_CRITICAL, "Error opening file '%s' for writing..\n\n", uvar.outputInjParams ); XLAL_ERROR ( XLAL_EIO ); } fprintf (fpInjParams, "%s", cfg.logString ); /* write search log comment */ if ( write_InjParams_to_fp ( fpInjParams, NULL, 0, uvar.outputMmunuX, numDetectors ) != XLAL_SUCCESS ) { /* write header-line comment */ XLAL_ERROR ( XLAL_EFUNC ); } } /* if outputInjParams */ multiAMBuffer_t XLAL_INIT_DECL(multiAMBuffer); /* prepare AM-buffer */ /* ----- prepare BSGL computation */ BSGLSetup *BSGLsetup = NULL; if ( uvar.computeBSGL ) { BOOLEAN useLogCorrection = TRUE; REAL4 *oLGX_p = NULL; REAL4 oLGX[PULSAR_MAX_DETECTORS]; if ( uvar.oLGX != NULL ) { XLAL_CHECK ( uvar.oLGX->length == numDetectors, XLAL_EINVAL, "Invalid input: length(oLGX) = %d differs from number of detectors (%d)'\n", uvar.oLGX->length, numDetectors ); XLAL_CHECK ( XLALParseLinePriors ( &oLGX[0], uvar.oLGX ) == XLAL_SUCCESS, XLAL_EFUNC ); oLGX_p = &oLGX[0]; } XLAL_CHECK ( ( BSGLsetup = XLALCreateBSGLSetup ( numDetectors, uvar.Fstar0, oLGX_p, useLogCorrection ) ) != NULL, XLAL_EFUNC ); } // if computeBSGL /* ----- main MC loop over numDraws trials ---------- */ INT4 i; for ( i=0; i < uvar.numDraws; i ++ ) { InjParams_t XLAL_INIT_DECL(injParamsDrawn); /* ----- generate signal random draws from ranges and generate Fstat atoms */ MultiFstatAtomVector *multiAtoms; multiAtoms = XLALSynthesizeTransientAtoms ( &injParamsDrawn, cfg.skypos, cfg.AmpPrior, cfg.transientInjectRange, cfg.multiDetStates, cfg.SignalOnly, &multiAMBuffer, cfg.rng, lineX, cfg.multiNoiseWeights ); XLAL_CHECK ( multiAtoms != NULL, XLAL_EFUNC ); /* ----- if requested, output signal injection parameters into file */ if ( fpInjParams && (write_InjParams_to_fp ( fpInjParams, &injParamsDrawn, uvar.dataStartGPS, uvar.outputMmunuX, numDetectors ) != XLAL_SUCCESS ) ) { XLAL_ERROR ( XLAL_EFUNC ); } /* if fpInjParams & failure*/ /* initialise BSGLComponents structure and allocate memory */ BSGLComponents XLAL_INIT_DECL(synthStats); /* struct containing multi-detector Fstat, single-detector Fstats, line-robust stat */ synthStats.numDetectors = numDetectors; /* compute F- and BSGListics from atoms */ UINT4 X; for ( X=0; X < numDetectors; X++ ) { synthStats.TwoFX[X] = XLALComputeFstatFromAtoms ( multiAtoms, X ); if ( xlalErrno != 0 ) { XLALPrintError ("\nError in function %s, line %d : Failed call to XLALComputeFstatFromAtoms().\n\n", __func__, __LINE__); XLAL_ERROR ( XLAL_EFUNC ); } } synthStats.TwoF = XLALComputeFstatFromAtoms ( multiAtoms, -1 ); if ( xlalErrno != 0 ) { XLALPrintError ("\nError in function %s, line %d : Failed call to XLALComputeFstatFromAtoms().\n\n", __func__, __LINE__); XLAL_ERROR ( XLAL_EFUNC ); } if ( uvar.computeBSGL ) { synthStats.log10BSGL = XLALComputeBSGL ( synthStats.TwoF, synthStats.TwoFX, BSGLsetup ); XLAL_CHECK ( xlalErrno == 0, XLAL_EFUNC, "XLALComputeBSGL() failed with xlalErrno = %d\n", xlalErrno ); } /* ----- if requested, output atoms-vector into file */ if ( uvar.outputAtoms ) { FILE *fpAtoms; char *fnameAtoms; UINT4 len = strlen ( uvar.outputAtoms ) + 20; if ( (fnameAtoms = XLALCalloc ( 1, len )) == NULL ) { XLALPrintError ("%s: failed to XLALCalloc ( 1, %d )\n", __func__, len ); XLAL_ERROR ( XLAL_EFUNC ); } sprintf ( fnameAtoms, "%s_%04d_of_%04d.dat", uvar.outputAtoms, i + 1, uvar.numDraws ); if ( ( fpAtoms = fopen ( fnameAtoms, "wb" )) == NULL ) { XLALPrintError ("%s: failed to open atoms-output file '%s' for writing.\n", __func__, fnameAtoms ); XLAL_ERROR ( XLAL_EFUNC ); } fprintf ( fpAtoms, "%s", cfg.logString ); /* output header info */ if ( write_MultiFstatAtoms_to_fp ( fpAtoms, multiAtoms ) != XLAL_SUCCESS ) { XLALPrintError ("%s: failed to write atoms to output file '%s'. xlalErrno = %d\n", __func__, fnameAtoms, xlalErrno ); XLAL_ERROR ( XLAL_EFUNC ); } XLALFree ( fnameAtoms ); fclose (fpAtoms); } /* if outputAtoms */ /* ----- if requested, output transient-cand statistics */ if ( fpStats && write_BSGL_candidate_to_fp ( fpStats, &synthStats, uvar.IFOs, &injParamsDrawn, uvar.computeBSGL ) != XLAL_SUCCESS ) { XLALPrintError ( "%s: write_transientCandidate_to_fp() failed.\n", __func__ ); XLAL_ERROR ( XLAL_EFUNC ); } /* ----- free Memory */ XLALDestroyMultiFstatAtomVector ( multiAtoms ); } /* for i < numDraws */ /* ----- close files ----- */ if ( fpStats ) fclose ( fpStats ); if ( fpInjParams ) fclose ( fpInjParams ); /* ----- free memory ---------- */ XLALDestroyMultiDetectorStateSeries ( cfg.multiDetStates ); XLALDestroyMultiNoiseWeights ( cfg.multiNoiseWeights ); XLALDestroyExpLUT(); XLALDestroyMultiAMCoeffs ( multiAMBuffer.multiAM ); /* ----- free amplitude prior pdfs ----- */ XLALDestroyPDF1D ( cfg.AmpPrior.pdf_h0Nat ); XLALDestroyPDF1D ( cfg.AmpPrior.pdf_cosi ); XLALDestroyPDF1D ( cfg.AmpPrior.pdf_psi ); XLALDestroyPDF1D ( cfg.AmpPrior.pdf_phi0 ); XLALFree ( BSGLsetup ); BSGLsetup = NULL; if ( cfg.logString ) { XLALFree ( cfg.logString ); } gsl_rng_free ( cfg.rng ); XLALDestroyUserVars(); /* did we forget anything ? (doesn't cover gsl-memory!) */ LALCheckMemoryLeaks(); return 0; } /* main() */
/** * MAIN function * Generates samples of B-stat and F-stat according to their pdfs for given signal-params. */ int main(int argc,char *argv[]) { UserInput_t XLAL_INIT_DECL(uvar); ConfigVariables XLAL_INIT_DECL(cfg); /**< various derived configuration settings */ vrbflg = 1; /* verbose error-messages */ LogSetLevel(lalDebugLevel); /* turn off default GSL error handler */ gsl_set_error_handler_off (); /* ----- register and read all user-variables ----- */ LogSetLevel(lalDebugLevel); if ( XLALInitUserVars( &uvar ) != XLAL_SUCCESS ) { LogPrintf ( LOG_CRITICAL, "%s: XLALInitUserVars() failed with errno=%d\n", __func__, xlalErrno ); return 1; } /* do ALL cmdline and cfgfile handling */ if ( XLALUserVarReadAllInput ( argc, argv ) != XLAL_SUCCESS ) { LogPrintf ( LOG_CRITICAL, "%s: XLALUserVarReadAllInput() failed with errno=%d\n", __func__, xlalErrno ); return 1; } if (uvar.help) /* if help was requested, we're done here */ return 0; if ( uvar.version ) { /* output verbose VCS version string if requested */ CHAR *vcs; if ( (vcs = XLALGetVersionString (lalDebugLevel)) == NULL ) { LogPrintf ( LOG_CRITICAL, "%s:XLALGetVersionString(%d) failed with errno=%d.\n", __func__, lalDebugLevel, xlalErrno ); return 1; } printf ( "%s\n", vcs ); XLALFree ( vcs ); return 0; } /* ---------- Initialize code-setup ---------- */ if ( XLALInitCode( &cfg, &uvar ) != XLAL_SUCCESS ) { LogPrintf (LOG_CRITICAL, "%s: XLALInitCode() failed with error = %d\n", __func__, xlalErrno ); XLAL_ERROR ( XLAL_EFUNC ); } /* ----- prepare stats output ----- */ FILE *fpTransientStats = NULL; if ( uvar.outputStats ) { if ( (fpTransientStats = fopen (uvar.outputStats, "wb")) == NULL) { LogPrintf (LOG_CRITICAL, "Error opening file '%s' for writing..\n\n", uvar.outputStats ); XLAL_ERROR ( XLAL_EIO ); } fprintf (fpTransientStats, "%s", cfg.logString ); /* write search log comment */ if ( write_transientCandidate_to_fp ( fpTransientStats, NULL ) != XLAL_SUCCESS ) { /* write header-line comment */ XLAL_ERROR ( XLAL_EFUNC ); } } /* if outputStats */ /* ----- prepare injection params output ----- */ FILE *fpInjParams = NULL; if ( uvar.outputInjParams ) { if ( (fpInjParams = fopen (uvar.outputInjParams, "wb")) == NULL) { LogPrintf (LOG_CRITICAL, "Error opening file '%s' for writing..\n\n", uvar.outputInjParams ); XLAL_ERROR ( XLAL_EIO ); } fprintf (fpInjParams, "%s", cfg.logString ); /* write search log comment */ if ( write_InjParams_to_fp ( fpInjParams, NULL, 0, 0, 0 ) != XLAL_SUCCESS ) { /* write header-line comment - options outputMmunuX and numDetectors not supported here, so pass defaults to deactivate them */ XLAL_ERROR ( XLAL_EFUNC ); } } /* if outputInjParams */ /* ----- main MC loop over numDraws trials ---------- */ multiAMBuffer_t XLAL_INIT_DECL(multiAMBuffer); /* prepare AM-buffer */ INT4 i; for ( i=0; i < uvar.numDraws; i ++ ) { InjParams_t XLAL_INIT_DECL(injParamsDrawn); /* ----- generate signal random draws from ranges and generate Fstat atoms */ MultiFstatAtomVector *multiAtoms; multiAtoms = XLALSynthesizeTransientAtoms ( &injParamsDrawn, cfg.skypos, cfg.AmpPrior, cfg.transientInjectRange, cfg.multiDetStates, cfg.SignalOnly, &multiAMBuffer, cfg.rng, -1, NULL ); // options lineX and noise_weights not supported here, so pass defaults to deactivate them if ( multiAtoms ==NULL ) { LogPrintf ( LOG_CRITICAL, "%s: XLALSynthesizeTransientAtoms() failed with xlalErrno = %d\n", __func__, xlalErrno ); XLAL_ERROR ( XLAL_EFUNC ); } /* ----- if requested, output signal injection parameters into file */ if ( fpInjParams && (write_InjParams_to_fp ( fpInjParams, &injParamsDrawn, uvar.dataStartGPS, 0, 0 ) ) != XLAL_SUCCESS ) { // options outputMmunuX and numDetectors not supported here, so pass defaults to deactivate them XLAL_ERROR ( XLAL_EFUNC ); } /* if fpInjParams & failure*/ /* ----- add meta-info on current transient-CW candidate */ transientCandidate_t XLAL_INIT_DECL(cand); cand.doppler.Alpha = multiAMBuffer.skypos.longitude; cand.doppler.Delta = multiAMBuffer.skypos.latitude; cand.windowRange = cfg.transientSearchRange; /* ----- if needed: compute transient-Bstat search statistic on these atoms */ if ( fpTransientStats || uvar.outputFstatMap || uvar.outputPosteriors ) { /* compute Fstat map F_mn over {t0, tau} */ if ( (cand.FstatMap = XLALComputeTransientFstatMap ( multiAtoms, cand.windowRange, uvar.useFReg)) == NULL ) { XLALPrintError ("%s: XLALComputeTransientFstatMap() failed with xlalErrno = %d.\n", __func__, xlalErrno ); XLAL_ERROR ( XLAL_EFUNC ); } } /* if we'll need the Fstat-map F_mn */ /* ----- if requested compute marginalized Bayes factor */ if ( fpTransientStats ) { cand.logBstat = XLALComputeTransientBstat ( cand.windowRange, cand.FstatMap ); UINT4 err = xlalErrno; if ( err ) { XLALPrintError ("%s: XLALComputeTransientBstat() failed with xlalErrno = %d\n", __func__, err ); XLAL_ERROR ( XLAL_EFUNC ); } if ( uvar.SignalOnly ) { cand.FstatMap->maxF += 2; cand.logBstat += 2; } } /* if Bstat requested */ /* ----- if requested, compute parameter posteriors for {t0, tau} */ pdf1D_t *pdf_t0 = NULL; pdf1D_t *pdf_tau = NULL; if ( fpTransientStats || uvar.outputPosteriors ) { if ( (pdf_t0 = XLALComputeTransientPosterior_t0 ( cand.windowRange, cand.FstatMap )) == NULL ) { XLALPrintError ("%s: failed to compute t0-posterior\n", __func__ ); XLAL_ERROR ( XLAL_EFUNC ); } if ( (pdf_tau = XLALComputeTransientPosterior_tau ( cand.windowRange, cand.FstatMap )) == NULL ) { XLALPrintError ("%s: failed to compute tau-posterior\n", __func__ ); XLAL_ERROR ( XLAL_EFUNC ); } /* get maximum-posterior estimate (MP) from the modes of these pdfs */ cand.t0_MP = XLALFindModeOfPDF1D ( pdf_t0 ); if ( xlalErrno ) { XLALPrintError ("%s: mode-estimation failed for pdf_t0. xlalErrno = %d\n", __func__, xlalErrno ); XLAL_ERROR ( XLAL_EFUNC ); } cand.tau_MP = XLALFindModeOfPDF1D ( pdf_tau ); if ( xlalErrno ) { XLALPrintError ("%s: mode-estimation failed for pdf_tau. xlalErrno = %d\n", __func__, xlalErrno ); XLAL_ERROR ( XLAL_EFUNC ); } } // if posteriors required /* ----- if requested, compute Ftotal over full data-span */ if ( uvar.computeFtotal ) { transientFstatMap_t *FtotalMap; /* prepare special window to cover all the data with one F-stat calculation == Ftotal */ transientWindowRange_t XLAL_INIT_DECL(winRangeAll); winRangeAll.type = TRANSIENT_NONE; BOOLEAN useFReg = false; if ( (FtotalMap = XLALComputeTransientFstatMap ( multiAtoms, winRangeAll, useFReg)) == NULL ) { XLALPrintError ("%s: XLALComputeTransientFstatMap() failed with xlalErrno = %d.\n", __func__, xlalErrno ); XLAL_ERROR ( XLAL_EFUNC ); } /* we only use twoFtotal = 2 * maxF from this single-Fstat calculation */ REAL8 twoFtotal = 2.0 * FtotalMap->maxF; if ( uvar.SignalOnly ) twoFtotal += 4; /* ugly hack: lacking a good container for twoFtotal, we borrow fkdot[3] for this here ;) [only used for paper-MCs] */ cand.doppler.fkdot[3] = twoFtotal; /* good riddance .. */ XLALDestroyTransientFstatMap ( FtotalMap ); } /* if computeFtotal */ /* ----- if requested, output atoms-vector into file */ if ( uvar.outputAtoms ) { FILE *fpAtoms; char *fnameAtoms; UINT4 len = strlen ( uvar.outputAtoms ) + 20; if ( (fnameAtoms = XLALCalloc ( 1, len )) == NULL ) { XLALPrintError ("%s: failed to XLALCalloc ( 1, %d )\n", __func__, len ); XLAL_ERROR ( XLAL_EFUNC ); } sprintf ( fnameAtoms, "%s_%04d_of_%04d.dat", uvar.outputAtoms, i + 1, uvar.numDraws ); if ( ( fpAtoms = fopen ( fnameAtoms, "wb" )) == NULL ) { XLALPrintError ("%s: failed to open atoms-output file '%s' for writing.\n", __func__, fnameAtoms ); XLAL_ERROR ( XLAL_EFUNC ); } fprintf ( fpAtoms, "%s", cfg.logString ); /* output header info */ if ( write_MultiFstatAtoms_to_fp ( fpAtoms, multiAtoms ) != XLAL_SUCCESS ) { XLALPrintError ("%s: failed to write atoms to output file '%s'. xlalErrno = %d\n", __func__, fnameAtoms, xlalErrno ); XLAL_ERROR ( XLAL_EFUNC ); } XLALFree ( fnameAtoms ); fclose (fpAtoms); } /* if outputAtoms */ /* ----- if requested, output Fstat-map over {t0, tau} */ if ( uvar.outputFstatMap ) { FILE *fpFstatMap; char *fnameFstatMap; UINT4 len = strlen ( uvar.outputFstatMap ) + 20; if ( (fnameFstatMap = XLALCalloc ( 1, len )) == NULL ) { XLALPrintError ("%s: failed to XLALCalloc ( 1, %d )\n", __func__, len ); XLAL_ERROR ( XLAL_EFUNC ); } sprintf ( fnameFstatMap, "%s_%04d_of_%04d.dat", uvar.outputFstatMap, i + 1, uvar.numDraws ); if ( ( fpFstatMap = fopen ( fnameFstatMap, "wb" )) == NULL ) { XLALPrintError ("%s: failed to open Fstat-map output file '%s' for writing.\n", __func__, fnameFstatMap ); XLAL_ERROR ( XLAL_EFUNC ); } fprintf ( fpFstatMap, "%s", cfg.logString ); /* output header info */ fprintf (fpFstatMap, "\nFstat_mn = \\\n" ); if ( XLALfprintfGSLmatrix ( fpFstatMap, "%.9g", cand.FstatMap->F_mn ) != XLAL_SUCCESS ) { XLALPrintError ("%s: XLALfprintfGSLmatrix() failed.\n", __func__ ); XLAL_ERROR ( XLAL_EFUNC ); } XLALFree ( fnameFstatMap ); fclose (fpFstatMap); } /* if outputFstatMap */ /* ----- if requested, output posterior pdfs on transient params {t0, tau} into a file */ if ( uvar.outputPosteriors ) { FILE *fpPosteriors; char *fnamePosteriors; UINT4 len = strlen ( uvar.outputPosteriors ) + 20; if ( (fnamePosteriors = XLALCalloc ( 1, len )) == NULL ) { XLALPrintError ("%s: failed to XLALCalloc ( 1, %d )\n", __func__, len ); XLAL_ERROR ( XLAL_EFUNC ); } sprintf ( fnamePosteriors, "%s_%04d_of_%04d.dat", uvar.outputPosteriors, i + 1, uvar.numDraws ); if ( ( fpPosteriors = fopen ( fnamePosteriors, "wb" )) == NULL ) { XLALPrintError ("%s: failed to open posteriors-output file '%s' for writing.\n", __func__, fnamePosteriors ); XLAL_ERROR ( XLAL_EFUNC ); } fprintf ( fpPosteriors, "%s", cfg.logString ); /* output header info */ /* write them to file, using pdf-method */ if ( XLALOutputPDF1D_to_fp ( fpPosteriors, pdf_t0, "pdf_t0" ) != XLAL_SUCCESS ) { XLALPrintError ("%s: failed to output t0-posterior to file '%s'.\n", __func__, fnamePosteriors ); XLAL_ERROR ( XLAL_EFUNC ); } if ( XLALOutputPDF1D_to_fp ( fpPosteriors, pdf_tau, "pdf_tau" ) != XLAL_SUCCESS ) { XLALPrintError ("%s: failed to output tau-posterior to file '%s'.\n", __func__, fnamePosteriors ); XLAL_ERROR ( XLAL_EFUNC ); } /* free mem, close file */ XLALFree ( fnamePosteriors ); fclose (fpPosteriors); } /* if outputPosteriors */ /* ----- if requested, output transient-cand statistics */ if ( fpTransientStats && write_transientCandidate_to_fp ( fpTransientStats, &cand ) != XLAL_SUCCESS ) { XLALPrintError ( "%s: write_transientCandidate_to_fp() failed.\n", __func__ ); XLAL_ERROR ( XLAL_EFUNC ); } /* ----- free Memory */ XLALDestroyTransientFstatMap ( cand.FstatMap ); XLALDestroyMultiFstatAtomVector ( multiAtoms ); XLALDestroyPDF1D ( pdf_t0 ); XLALDestroyPDF1D ( pdf_tau ); } /* for i < numDraws */ /* ----- close files ----- */ if ( fpTransientStats) fclose ( fpTransientStats ); if ( fpInjParams ) fclose ( fpInjParams ); /* ----- free memory ---------- */ XLALDestroyMultiDetectorStateSeries ( cfg.multiDetStates ); XLALDestroyMultiAMCoeffs ( multiAMBuffer.multiAM ); XLALDestroyExpLUT(); /* ----- free amplitude prior pdfs ----- */ XLALDestroyPDF1D ( cfg.AmpPrior.pdf_h0Nat ); XLALDestroyPDF1D ( cfg.AmpPrior.pdf_cosi ); XLALDestroyPDF1D ( cfg.AmpPrior.pdf_psi ); XLALDestroyPDF1D ( cfg.AmpPrior.pdf_phi0 ); if ( cfg.logString ) XLALFree ( cfg.logString ); gsl_rng_free ( cfg.rng ); XLALDestroyUserVars(); /* did we forget anything ? (doesn't cover gsl-memory!) */ LALCheckMemoryLeaks(); return 0; } /* main() */
/** * Compute transient-CW posterior (normalized) on timescale tau, using given type and parameters * of transient window range. * * NOTE: the returned pdf has a number of sample-points N_tauRange given by the size * of the input matrix FstatMap (namely N_tauRange = tauBand / dtau) * */ pdf1D_t * XLALComputeTransientPosterior_tau ( transientWindowRange_t windowRange, /**< [in] type and parameters specifying transient window range */ const transientFstatMap_t *FstatMap /**< [in] pre-computed transient-Fstat map F_mn over {t0, tau} ranges */ ) { /* ----- check input consistency */ if ( !FstatMap || !FstatMap->F_mn ) { XLALPrintError ("%s: invalid NULL input 'FstatMap' or 'FstatMap->F_mn'\n", __func__ ); XLAL_ERROR_NULL ( XLAL_EINVAL ); } if ( windowRange.type >= TRANSIENT_LAST ) { XLALPrintError ("%s: unknown window-type (%d) passes as input. Allowed are [0,%d].\n", __func__, windowRange.type, TRANSIENT_LAST-1); XLAL_ERROR_NULL ( XLAL_EINVAL ); } /* ----- step through F_mn array subtract maxF and sum e^{F_mn - maxF}*/ /* * It is numerically more robust to marginalize over e^(F_mn - Fmax), which at worst can underflow, while * e^F_mn can overflow (for F>~700). The constant offset e^Fmax is irrelevant for posteriors (normalization constant). */ UINT4 N_t0Range = FstatMap->F_mn->size1; UINT4 N_tauRange = FstatMap->F_mn->size2; REAL8 tau0 = windowRange.tau; REAL8 tau1 = tau0 + windowRange.tauBand; pdf1D_t *ret; /* ----- handle special cases: 1) point-like support, 2) uniform pdf-value over 1 bin ----- */ if ( N_tauRange == 1 && (windowRange.tauBand == 0) ) { if ( (ret = XLALCreateSingularPDF1D ( tau0 )) == NULL ) { XLALPrintError ("%s: failed to create singular pdf for tau0 = %g\n", __func__, tau0 ); XLAL_ERROR_NULL ( XLAL_EFUNC ); } return ret; } /* if singular pdf in tau */ if ( (N_tauRange == 1) && (windowRange.tauBand > 0) ) { if ( (ret = XLALCreateUniformPDF1D ( tau0, tau1 )) == NULL ) { XLALPrintError ( "%s: failed to created unform pdf over [%g, %g]\n", __func__, tau0, tau1 ); XLAL_ERROR_NULL ( XLAL_EFUNC ); } return ret; } /* if uniform pdf over small band tauBand */ /* ----- general N>1 point pdf case ----- */ if ( ( ret = XLALCreateDiscretePDF1D ( tau0, tau1, N_tauRange )) == NULL ) { XLALPrintError ("%s: XLALCreateDiscretePDF1D() failed with xlalErrno = %d\n", __func__, xlalErrno ); XLAL_ERROR_NULL ( XLAL_ENOMEM ); } UINT4 m, n; for ( n=0; n < N_tauRange; n ++ ) /* loop over timescales tau */ { REAL8 sum_eF = 0; for ( m=0; m < N_t0Range; m ++ ) /* loop over start-times t0 */ { REAL8 DeltaF = FstatMap->maxF - gsl_matrix_get ( FstatMap->F_mn, m, n ); // always >= 0, exactly ==0 at {m,n}_max //sum_eB += exp ( - DeltaF ); sum_eF += XLALFastNegExp ( DeltaF ); } /* for m < N_t0Range */ ret->probDens->data[n] = sum_eF; } /* for n < N_tauRange */ /* free mem */ XLALDestroyExpLUT(); /* normalize this PDF */ if ( XLALNormalizePDF1D ( ret ) != XLAL_SUCCESS ) { XLALPrintError ("%s: failed to normalize posterior pdf ..\n", __func__ ); XLAL_ERROR_NULL ( XLAL_EFUNC ); } /* ----- return ----- */ return ret; } /* XLALComputeTransientPosterior_tau() */