REAL8FrequencySeries *XLALFrStreamInputREAL8FrequencySeries(LALFrStream * stream, const char *chname, const LIGOTimeGPS * epoch) { REAL8FrequencySeries *series; LALTYPECODE typecode; if (XLALFrStreamSeek(stream, epoch)) XLAL_ERROR_NULL(XLAL_EFUNC); typecode = XLALFrFileQueryChanType(stream->file, chname, stream->pos); switch (typecode) { case LAL_S_TYPE_CODE: INPUTFS(series, REAL8, REAL4, S2S, stream, chname, epoch); break; case LAL_D_TYPE_CODE: series = XLALFrStreamReadREAL8FrequencySeries(stream, chname, epoch); if (!series) XLAL_ERROR_NULL(XLAL_EFUNC); break; case LAL_C_TYPE_CODE: case LAL_Z_TYPE_CODE: XLAL_PRINT_ERROR("Cannot convert complex type to float type"); default: XLAL_ERROR_NULL(XLAL_ETYPE); } return series; }
/** * @brief Creates an equation of state structure from tabulated equation * of state data of a known name. * @details A known, installed, named tabulated equation of state data file is * read and the used to create the equation of state structure. Presently * the known equations of state are: * - AP4 * - FPS * - SLY4 * @param[in] name The name of the equation of state. * @return A pointer to neutron star equation of state structure. */ LALSimNeutronStarEOS *XLALSimNeutronStarEOSByName(const char *name) { static const char fname_base[] = "LALSimNeutronStarEOS_"; static const char fname_extn[] = ".dat"; static const char *eos_names[] = { "FPS", "SLY4", "AP4" }; size_t n = XLAL_NUM_ELEM(eos_names); size_t i; char fname[FILENAME_MAX]; for (i = 0; i < n; ++i) if (mystrcasecmp(name, eos_names[i]) == 0) { LALSimNeutronStarEOS *eos; snprintf(fname, sizeof(fname), "%s%s%s", fname_base, eos_names[i], fname_extn); eos = XLALSimNeutronStarEOSFromFile(fname); if (!eos) XLAL_ERROR_NULL(XLAL_EFUNC); snprintf(eos->name, sizeof(eos->name), "%s", eos_names[i]); return eos; } XLAL_PRINT_ERROR("Unrecognized EOS name %s...", name); XLALPrintError("\tKnown EOS names are: %s", eos_names[0]); for (i = 1; i < n; ++i) XLALPrintError(", %s", eos_names[i]); XLALPrintError("\n"); XLAL_ERROR_NULL(XLAL_ENAME); }
REAL8TimeSeries *XLALFrStreamInputREAL8TimeSeries(LALFrStream * stream, const char *chname, const LIGOTimeGPS * start, double duration, size_t lengthlimit) { REAL8TimeSeries *series; LALTYPECODE typecode; if (XLALFrStreamSeek(stream, start)) XLAL_ERROR_NULL(XLAL_EFUNC); typecode = XLALFrFileQueryChanType(stream->file, chname, stream->pos); switch (typecode) { case LAL_I2_TYPE_CODE: INPUTTS(series, REAL8, INT2, S2S, stream, chname, start, duration, lengthlimit); break; case LAL_I4_TYPE_CODE: INPUTTS(series, REAL8, INT4, S2S, stream, chname, start, duration, lengthlimit); break; case LAL_I8_TYPE_CODE: INPUTTS(series, REAL8, INT8, S2S, stream, chname, start, duration, lengthlimit); break; case LAL_U2_TYPE_CODE: INPUTTS(series, REAL8, UINT2, S2S, stream, chname, start, duration, lengthlimit); break; case LAL_U4_TYPE_CODE: INPUTTS(series, REAL8, UINT4, S2S, stream, chname, start, duration, lengthlimit); break; case LAL_U8_TYPE_CODE: INPUTTS(series, REAL8, UINT8, S2S, stream, chname, start, duration, lengthlimit); break; case LAL_S_TYPE_CODE: INPUTTS(series, REAL8, REAL4, S2S, stream, chname, start, duration, lengthlimit); break; case LAL_D_TYPE_CODE: series = XLALFrStreamReadREAL8TimeSeries(stream, chname, start, duration, lengthlimit); if (!series) XLAL_ERROR_NULL(XLAL_EFUNC); break; case LAL_C_TYPE_CODE: case LAL_Z_TYPE_CODE: XLAL_PRINT_ERROR("Cannot convert complex type to float type"); default: XLAL_ERROR_NULL(XLAL_ETYPE); } return series; }
/** * Parse string-vectors (typically input by user) of N detector noise-floors \f$\sqrt{S_X}\f$ * for detectors \f$X=1\ldots N\f$, where here we assume equal number of SFTs per detector * such that \f$S^{-1} = \frac{1}{N}\sum_{X=0}^{N-1} S_X^{-1}\f$. * * The detectors corresponding to each noise-floor may be a subset of the input string-vectors, * e.g. if parsing noise-floors for a segment where SFTs from some detectors are missing. * The vector \p */ int XLALParseMultiNoiseFloorMapped ( MultiNoiseFloor *multiNoiseFloor, /**< [out] parsed multi-IFO noise floor info */ const LALStringVector *multiNoiseFloorDetNames, /**< [in] detector names for entries in \p multiNoiseFloor */ const LALStringVector *sqrtSX, /**< [in] string-list of \f$\sqrt{S_X}\f$ for detectors \f$X\f$ */ const LALStringVector *sqrtSXDetNames /**< [in] detector names for entries in \p sqrtSX */ ) { XLAL_CHECK ( multiNoiseFloor != NULL, XLAL_EINVAL ); XLAL_CHECK ( multiNoiseFloorDetNames != NULL, XLAL_EINVAL ); XLAL_CHECK ( sqrtSX != NULL, XLAL_EINVAL ); XLAL_CHECK ( sqrtSXDetNames != NULL, XLAL_EINVAL ); XLAL_CHECK ( multiNoiseFloorDetNames->length <= sqrtSXDetNames->length, XLAL_EINVAL ); XLAL_CHECK ( sqrtSX->length == sqrtSXDetNames->length, XLAL_EINVAL ); /* parse input strings */ REAL8 sqrtSn[PULSAR_MAX_DETECTORS]; for ( UINT4 Y = 0; Y < sqrtSX->length; Y ++ ) { XLAL_CHECK ( XLALParseStringValueAsREAL8 ( &sqrtSn[Y], sqrtSX->data[Y] ) == XLAL_SUCCESS, XLAL_EFUNC ); XLAL_CHECK ( sqrtSn[Y] >= 0, XLAL_EDOM ); } /* initialize empty return struct */ multiNoiseFloor->length = multiNoiseFloorDetNames->length; /* fill multiNoiseFloor with correctly mapped values */ for ( UINT4 X = 0; X < multiNoiseFloor->length; X ++ ) { const INT4 Y = XLALFindStringInVector( multiNoiseFloorDetNames->data[X], sqrtSXDetNames ); if ( Y < 0 ) { char *sqrtSXDet = XLALConcatStringVector( sqrtSXDetNames, "," ); XLAL_PRINT_ERROR ( "Noise-floor detector '%s' not found in list of sqrtSX detectors '%s'", multiNoiseFloorDetNames->data[X], sqrtSXDet ); XLALFree ( sqrtSXDet ); XLAL_ERROR ( XLAL_EINVAL ); } multiNoiseFloor->sqrtSn[X] = sqrtSn[Y]; } return XLAL_SUCCESS; } /* XLALParseMultiNoiseFloorMapped() */
void LALInferenceTemplateXLALSimInspiralChooseWaveform(LALInferenceModel *model) /*************************************************************************************************************************/ /* Wrapper for LALSimulation waveforms: */ /* XLALSimInspiralChooseFDWaveform() and XLALSimInspiralChooseTDWaveform(). */ /* */ /* model->params parameters are: */ /* - "name" description; type OPTIONAL (default value) */ /* */ /* MODEL PARAMETERS */ /* - "LAL_APPROXIMANT Approximant; Approximant */ /* - "LAL_PNORDER" Phase PN order; INT4 */ /* - "LAL_AMPORDER" Amplitude PN order; INT4 OPTIONAL (-1) */ /* - "spinO" Spin order; LALSimInspiralSpinOrder OPTIONAL (LAL_SIM_INSPIRAL_SPIN_ORDER_DEFAULT) */ /* - "tideO" Tidal order; LALSimInspiralTidalOrder OPTIONAL (LAL_SIM_INSPIRAL_TIDAL_ORDER_DEFAULT) */ /* - "f_ref" frequency at which the (frequency dependent) parameters are defined; REAL8 OPTIONAL (0.0) */ /* - "fLow" lower frequency bound; REAL8 OPTIONAL (model->fLow) */ /* */ /* MASS PARAMETERS; either: */ /* - "mass1" mass of object 1 in solar mass; REAL8 */ /* - "mass2" mass of object 1 in solar mass; REAL8 */ /* OR */ /* - "chirpmass" chirpmass in solar mass; REAL8 */ /* - "q" asymmetric mass ratio m2/m1, 0<q<1; REAL8 */ /* OR */ /* - "chirpmass" chirpmass in solar mass; REAL8 */ /* - "eta" symmetric mass ratio (m1*m2)/(m1+m2)^2; REAL8 */ /* */ /* ORIENTATION AND SPIN PARAMETERS */ /* - "phi0" reference phase as per LALSimulation convention; REAL8 */ /* - "distance" distance in Mpc */ /* - "costheta_jn"); cos of zenith angle between J and N in radians; REAL8 */ /* - "phi_jl"); azimuthal angle of L_N on its cone about J radians; REAL8 */ /* - "tilt_spin1"); zenith angle between S1 and LNhat in radians; REAL8 */ /* - "tilt_spin2"); zenith angle between S2 and LNhat in radians; REAL8 */ /* - "phi12"); difference in azimuthal angle between S1, S2 in radians; REAL8 */ /* - "a_spin1" magnitude of spin 1 in general configuration, -1<a_spin1<1; REAL8 OPTIONAL (0.0) */ /* - "a_spin2" magnitude of spin 2 in general configuration, -1<a_spin1<1; REAL8 OPTIONAL (0.0) */ /* */ /* OTHER PARAMETERS */ /* - "lambda1" tidal parameter of object 1; REAL8 OPTIONAL (0.0) */ /* - "lambda2" tidal parameter of object 1; REAL8 OPTIONAL (0.0) */ /* */ /* - "time" used as an OUTPUT only; REAL8 */ /* */ /* */ /* model needs to also contain: */ /* - model->fLow Unless - "fLow" OPTIONAL */ /* - model->deltaT */ /* - if model->domain == LAL_SIM_DOMAIN_FREQUENCY */ /* - model->deltaF */ /* - model->freqhCross */ /* - model->freqhPlus */ /* - else */ /* - model->timehPlus */ /* - model->timehCross */ /*************************************************************************************************************************/ { Approximant approximant = (Approximant) 0; INT4 order=-1; INT4 amporder; static int sizeWarning = 0; int ret=0; INT4 errnum=0; REAL8TimeSeries *hplus=NULL; /**< +-polarization waveform [returned] */ REAL8TimeSeries *hcross=NULL; /**< x-polarization waveform [returned] */ COMPLEX16FrequencySeries *hptilde=NULL, *hctilde=NULL; REAL8 mc; REAL8 phi0, deltaT, m1, m2, f_low, f_start, distance, inclination; REAL8 *m1_p,*m2_p; REAL8 deltaF, f_max; /* Sampling rate for time domain models */ deltaT = model->deltaT; if (LALInferenceCheckVariable(model->params, "LAL_APPROXIMANT")) approximant = *(Approximant*) LALInferenceGetVariable(model->params, "LAL_APPROXIMANT"); else { XLALPrintError(" ERROR in templateLALGenerateInspiral(): (INT4) \"LAL_APPROXIMANT\" parameter not provided!\n"); XLAL_ERROR_VOID(XLAL_EDATA); } if (LALInferenceCheckVariable(model->params, "LAL_PNORDER")) order = *(INT4*) LALInferenceGetVariable(model->params, "LAL_PNORDER"); else { XLALPrintError(" ERROR in templateLALGenerateInspiral(): (INT4) \"LAL_PNORDER\" parameter not provided!\n"); XLAL_ERROR_VOID(XLAL_EDATA); } /* Explicitly set the default amplitude order if one is not specified. * This serves two purposes: * 1) The default behavior of the code won't change unexpectedly due to changes in LALSimulation. * 2) We need to know the amplitude order in order to set the starting frequency of the waveform properly. */ if (LALInferenceCheckVariable(model->params, "LAL_AMPORDER")) amporder = *(INT4*) LALInferenceGetVariable(model->params, "LAL_AMPORDER"); else amporder = -1; REAL8 f_ref = 100.0; if (LALInferenceCheckVariable(model->params, "f_ref")) f_ref = *(REAL8 *)LALInferenceGetVariable(model->params, "f_ref"); REAL8 fTemp = f_ref; if(LALInferenceCheckVariable(model->params,"chirpmass")) { mc = *(REAL8*) LALInferenceGetVariable(model->params, "chirpmass"); if (LALInferenceCheckVariable(model->params,"q")) { REAL8 q = *(REAL8 *)LALInferenceGetVariable(model->params,"q"); q2masses(mc, q, &m1, &m2); } else { REAL8 eta = *(REAL8*) LALInferenceGetVariable(model->params, "eta"); mc2masses(mc, eta, &m1, &m2); } } else if((m1_p=(REAL8 *)LALInferenceGetVariable(model->params, "mass1")) && (m2_p=(REAL8 *)LALInferenceGetVariable(model->params, "mass2"))) { m1=*m1_p; m2=*m2_p; } else { fprintf(stderr,"No mass parameters found!"); exit(0); } distance = LALInferenceGetREAL8Variable(model->params,"distance")* LAL_PC_SI * 1.0e6; /* distance (1 Mpc) in units of metres */ phi0 = LALInferenceGetREAL8Variable(model->params, "phase"); /* START phase as per lalsimulation convention, radians*/ /* Zenith angle between J and N in radians. Also known as inclination angle when spins are aligned */ REAL8 thetaJN = acos(LALInferenceGetREAL8Variable(model->params, "costheta_jn")); /* zenith angle between J and N in radians */ /* Check if fLow is a model parameter, otherwise use data structure definition */ if(LALInferenceCheckVariable(model->params, "flow")) f_low = *(REAL8*) LALInferenceGetVariable(model->params, "flow"); else f_low = model->fLow; f_start = fLow2fStart(f_low, amporder, approximant); f_max = 0.0; /* for freq domain waveforms this will stop at ISCO. Previously found using model->fHigh causes NaNs in waveform (see redmine issue #750)*/ /* ==== SPINS ==== */ /* We will default to spinless signal and then add in the spin components if required */ /* If there are non-aligned spins, we must convert between the System Frame coordinates * and the cartestian coordinates */ /* The cartesian spin coordinates (default 0), as passed to LALSimulation */ REAL8 spin1x = 0.0; REAL8 spin1y = 0.0; REAL8 spin1z = 0.0; REAL8 spin2x = 0.0; REAL8 spin2y = 0.0; REAL8 spin2z = 0.0; /* System frame coordinates as used for jump proposals */ REAL8 a_spin1 = 0.0; /* Magnitude of spin1 */ REAL8 a_spin2 = 0.0; /* Magnitude of spin2 */ REAL8 phiJL = 0.0; /* azimuthal angle of L_N on its cone about J radians */ REAL8 tilt1 = 0.0; /* zenith angle between S1 and LNhat in radians */ REAL8 tilt2 = 0.0; /* zenith angle between S2 and LNhat in radians */ REAL8 phi12 = 0.0; /* difference in azimuthal angle btwn S1, S2 in radians */ /* Now check if we have spin amplitudes */ if(LALInferenceCheckVariable(model->params, "a_spin1")) a_spin1 = *(REAL8*) LALInferenceGetVariable(model->params, "a_spin1"); if(LALInferenceCheckVariable(model->params, "a_spin2")) a_spin2 = *(REAL8*) LALInferenceGetVariable(model->params, "a_spin2"); /* Check if we have spin angles too */ if(LALInferenceCheckVariable(model->params, "phi_jl")) phiJL = LALInferenceGetREAL8Variable(model->params, "phi_jl"); if(LALInferenceCheckVariable(model->params, "tilt_spin1")) tilt1 = LALInferenceGetREAL8Variable(model->params, "tilt_spin1"); if(LALInferenceCheckVariable(model->params, "tilt_spin2")) tilt2 = LALInferenceGetREAL8Variable(model->params, "tilt_spin2"); if(LALInferenceCheckVariable(model->params, "phi12")) phi12 = LALInferenceGetREAL8Variable(model->params, "phi12"); /* If we have tilt angles zero, then the spins are aligned and we just set the z component */ /* However, if the waveform supports precession then we still need to get the right coordinate components */ SpinSupport spin_support=XLALSimInspiralGetSpinSupportFromApproximant(approximant); if(tilt1==0.0 && tilt2==0.0 && (spin_support==LAL_SIM_INSPIRAL_SPINLESS || spin_support==LAL_SIM_INSPIRAL_ALIGNEDSPIN)) { spin1z=a_spin1; spin2z=a_spin2; inclination = thetaJN; /* Inclination angle is just thetaJN */ } else { /* Template is not aligned-spin only. */ /* Set all the other spin components according to the angles we received above */ /* The transformation function doesn't know fLow, so f_ref==0 isn't interpretted as a request to use the starting frequency for reference. */ if(fTemp==0.0) fTemp = f_start; XLAL_TRY(ret=XLALSimInspiralTransformPrecessingNewInitialConditions( &inclination, &spin1x, &spin1y, &spin1z, &spin2x, &spin2y, &spin2z, thetaJN, phiJL, tilt1, tilt2, phi12, a_spin1, a_spin2, m1*LAL_MSUN_SI, m2*LAL_MSUN_SI, fTemp), errnum); if (ret == XLAL_FAILURE) { XLALPrintError(" ERROR in XLALSimInspiralTransformPrecessingNewInitialConditions(): error converting angles. errnum=%d\n",errnum ); return; } } /* ==== TIDAL PARAMETERS ==== */ REAL8 lambda1 = 0.; if(LALInferenceCheckVariable(model->params, "lambda1")) lambda1 = *(REAL8*) LALInferenceGetVariable(model->params, "lambda1"); REAL8 lambda2 = 0.; if(LALInferenceCheckVariable(model->params, "lambda2")) lambda2 = *(REAL8*) LALInferenceGetVariable(model->params, "lambda2"); REAL8 lambdaT = 0.; REAL8 dLambdaT = 0.; REAL8 sym_mass_ratio_eta = 0.; if(LALInferenceCheckVariable(model->params, "lambdaT")&&LALInferenceCheckVariable(model->params, "dLambdaT")){ lambdaT = *(REAL8*) LALInferenceGetVariable(model->params, "lambdaT"); dLambdaT = *(REAL8*) LALInferenceGetVariable(model->params, "dLambdaT"); sym_mass_ratio_eta = m1*m2/((m1+m2)*(m1+m2)); LALInferenceLambdaTsEta2Lambdas(lambdaT,dLambdaT,sym_mass_ratio_eta,&lambda1,&lambda2); } /* Only use GR templates */ LALSimInspiralTestGRParam *nonGRparams = NULL; /* ==== Call the waveform generator ==== */ if(model->domain == LAL_SIM_DOMAIN_FREQUENCY) { deltaF = model->deltaF; XLAL_TRY(ret=XLALSimInspiralChooseFDWaveformFromCache(&hptilde, &hctilde, phi0, deltaF, m1*LAL_MSUN_SI, m2*LAL_MSUN_SI, spin1x, spin1y, spin1z, spin2x, spin2y, spin2z, f_start, f_max, f_ref, distance, inclination,lambda1, lambda2, model->waveFlags, nonGRparams, amporder, order, approximant,model->waveformCache, NULL), errnum); /* if the waveform failed to generate, fill the buffer with zeros * so that the previous waveform is not left there */ if (ret != XLAL_SUCCESS || hptilde == NULL || hctilde == NULL) { XLALPrintError(" ERROR in XLALSimInspiralChooseWaveformFromCache(): error generating waveform. errnum=%d\n",errnum ); memset(model->freqhPlus->data->data,0,sizeof(model->freqhPlus->data->data[0])*model->freqhPlus->data->length); memset(model->freqhCross->data->data,0,sizeof(model->freqhCross->data->data[0])*model->freqhCross->data->length); if ( hptilde ) XLALDestroyCOMPLEX16FrequencySeries(hptilde); if ( hctilde ) XLALDestroyCOMPLEX16FrequencySeries(hctilde); XLALSimInspiralDestroyTestGRParam(nonGRparams); XLAL_ERROR_VOID(XLAL_FAILURE); } if (hptilde==NULL || hptilde->data==NULL || hptilde->data->data==NULL ) { XLALPrintError(" ERROR in LALInferenceTemplateXLALSimInspiralChooseWaveform(): encountered unallocated 'hptilde'.\n"); XLAL_ERROR_VOID(XLAL_EFAULT); } if (hctilde==NULL || hctilde->data==NULL || hctilde->data->data==NULL ) { XLALPrintError(" ERROR in LALInferenceTemplateXLALSimInspiralChooseWaveform(): encountered unallocated 'hctilde'.\n"); XLAL_ERROR_VOID(XLAL_EFAULT); } INT4 rem=0; UINT4 size=hptilde->data->length; if(size>model->freqhPlus->data->length) size=model->freqhPlus->data->length; memcpy(model->freqhPlus->data->data,hptilde->data->data,sizeof(hptilde->data->data[0])*size); if( (rem=(model->freqhPlus->data->length - size)) > 0) memset(&(model->freqhPlus->data->data[size]),0, rem*sizeof(hptilde->data->data[0]) ); size=hctilde->data->length; if(size>model->freqhCross->data->length) size=model->freqhCross->data->length; memcpy(model->freqhCross->data->data,hctilde->data->data,sizeof(hctilde->data->data[0])*size); if( (rem=(model->freqhCross->data->length - size)) > 0) memset(&(model->freqhCross->data->data[size]),0, rem*sizeof(hctilde->data->data[0]) ); /* Destroy the nonGr params */ XLALSimInspiralDestroyTestGRParam(nonGRparams); REAL8 instant = model->freqhPlus->epoch.gpsSeconds + 1e-9*model->freqhPlus->epoch.gpsNanoSeconds; LALInferenceSetVariable(model->params, "time", &instant); } else { XLAL_TRY(ret=XLALSimInspiralChooseTDWaveformFromCache(&hplus, &hcross, phi0, deltaT, m1*LAL_MSUN_SI, m2*LAL_MSUN_SI, spin1x, spin1y, spin1z, spin2x, spin2y, spin2z, f_start, f_ref, distance, inclination, lambda1, lambda2, model->waveFlags, nonGRparams, amporder, order, approximant,model->waveformCache), errnum); XLALSimInspiralDestroyTestGRParam(nonGRparams); if (ret == XLAL_FAILURE || hplus == NULL || hcross == NULL) { XLALPrintError(" ERROR in XLALSimInspiralChooseWaveformFromCache(): error generating waveform. errnum=%d\n",errnum ); memset(model->timehPlus->data->data,0,sizeof(model->timehPlus->data->data[0]) * model->timehPlus->data->length); memset(model->timehCross->data->data,0,sizeof(model->timehCross->data->data[0]) * model->timehCross->data->length); if ( hplus ) XLALDestroyREAL8TimeSeries(hplus); if ( hcross ) XLALDestroyREAL8TimeSeries(hcross); XLALSimInspiralDestroyTestGRParam(nonGRparams); XLAL_ERROR_VOID(XLAL_FAILURE); } /* The following complicated mess is a result of the following considerations: 1) The discrete time samples of the template and the timeModel buffers will not, in general line up. 2) The likelihood function will timeshift the template in the frequency domain to align it properly with the desired tc in each detector (these are different because the detectors receive the signal at different times). Because this timeshifting is done in the frequency domain, the effective time-domain template is periodic. We want to avoid the possibility of non-zero template samples wrapping around from the start/end of the buffer, since real templates are not periodic! 3) If the template apporaches the ends of the timeModel buffer, then it should be tapered in the same way as the timeData (currently 0.4 seconds, hard-coded! Tukey window; see LALInferenceReadData.c, near line 233) so that template and signal in the data match. However, as an optimization, we perform only one tapering and FFT-ing in the likelihood function; subsequent timeshifts for the different detectors will cause the tapered regions of the template and data to become mis-aligned. The algorthim we use is the following: 1) Inject the template to align with the nearest sample in the timeModel buffer to the desired geocent_end time. 2) Check whether either the start or the end of the template overlaps the tapered region, plus a safety buffer corresponding to a conservative estimate of the largest geocenter <--> detector timeshift. a) If there is no overlap at the start or end of the buffer, we're done. b) If there is an overlap, issue one warning per process (which can be disabled by setting the LAL debug level) about a too-short segment length, and return. */ size_t waveLength = hplus->data->length; size_t bufLength = model->timehPlus->data->length; /* 2*Rearth/(c*deltaT)---2 is safety factor---is the maximum time shift for any earth-based detector. */ size_t maxShift = (size_t)lround(4.255e-2/deltaT); /* Taper 0.4 seconds at start and end (hard-coded! in LALInferenceReadData.c, around line 233). */ size_t taperLength = (size_t)lround(0.4/deltaT); /* Within unsafeLength of ends of buffer, possible danger of wrapping and/or tapering interactions. */ size_t unsafeLength = taperLength + maxShift; REAL8 desiredTc = *(REAL8 *)LALInferenceGetVariable(model->params, "time"); REAL8 tStart = XLALGPSGetREAL8(&(model->timehPlus->epoch)); REAL8 tEnd = tStart + deltaT * model->timehPlus->data->length; if (desiredTc < tStart || desiredTc > tEnd) { XLALDestroyREAL8TimeSeries(hplus); XLALDestroyREAL8TimeSeries(hcross); XLAL_PRINT_ERROR("desired tc (%.4f) outside data buffer\n", desiredTc); XLAL_ERROR_VOID(XLAL_EDOM); } /* The nearest sample in model buffer to the desired tc. */ size_t tcSample = (size_t)lround((desiredTc - XLALGPSGetREAL8(&(model->timehPlus->epoch)))/deltaT); /* The acutal coalescence time that corresponds to the buffer sample on which the waveform's tC lands. */ REAL8 injTc = XLALGPSGetREAL8(&(model->timehPlus->epoch)) + tcSample*deltaT; /* The sample at which the waveform reaches tc. */ size_t waveTcSample = (size_t)lround(-XLALGPSGetREAL8(&(hplus->epoch))/deltaT); /* 1 + (number of samples post-tc in waveform) */ size_t wavePostTc = waveLength - waveTcSample; size_t bufStartIndex = (tcSample >= waveTcSample ? tcSample - waveTcSample : 0); size_t bufEndIndex = (wavePostTc + tcSample <= bufLength ? wavePostTc + tcSample : bufLength); size_t bufWaveLength = bufEndIndex - bufStartIndex; size_t waveStartIndex = (tcSample >= waveTcSample ? 0 : waveTcSample - tcSample); if (bufStartIndex < unsafeLength || (bufLength - bufEndIndex) <= unsafeLength) { /* The waveform could be timeshifted into a region where it will be tapered improperly, or even wrap around from the periodic timeshift. Issue warning. */ if (!sizeWarning) { fprintf(stderr, "WARNING: Generated template is too long to guarantee that it will not\n"); fprintf(stderr, "WARNING: (a) lie in a tapered region of the time-domain buffer\n"); fprintf(stderr, "WARNING: (b) wrap periodically when timeshifted in likelihood computation\n"); fprintf(stderr, "WARNING: Either of these may cause differences between the template and the\n"); fprintf(stderr, "WARNING: correct GW waveform in each detector.\n"); fprintf(stderr, "WARNING: Parameter estimation will continue, but you should consider\n"); fprintf(stderr, "WARNING: increasing the data segment length (using the --seglen) option.\n"); sizeWarning = 1; } } /* Clear model buffers */ memset(model->timehPlus->data->data, 0, sizeof(REAL8)*model->timehPlus->data->length); memset(model->timehCross->data->data, 0, sizeof(REAL8)*model->timehCross->data->length); /* Inject */ memcpy(model->timehPlus->data->data + bufStartIndex, hplus->data->data + waveStartIndex, bufWaveLength*sizeof(REAL8)); memcpy(model->timehCross->data->data + bufStartIndex, hcross->data->data + waveStartIndex, bufWaveLength*sizeof(REAL8)); LALInferenceSetVariable(model->params, "time", &injTc); } if ( hplus ) XLALDestroyREAL8TimeSeries(hplus); if ( hcross ) XLALDestroyREAL8TimeSeries(hcross); if ( hptilde ) XLALDestroyCOMPLEX16FrequencySeries(hptilde); if ( hctilde ) XLALDestroyCOMPLEX16FrequencySeries(hctilde); return; }
/** * \brief Gzip the nested sample files * * This function gzips the output nested sample files. It will also strip unneccesary parameters * from the header "_params.txt" file. * * \param runState [in] The analysis information structure */ void gzip_output( LALInferenceRunState *runState ){ /* Single thread here */ LALInferenceThreadState *threadState = runState->threads[0]; /* Open original output output file */ CHAR *outfile = NULL; ProcessParamsTable *ppt1 = LALInferenceGetProcParamVal( runState->commandLine, "--outfile" ); if( ppt1 ){ CHAR outfilepars[256] = "", outfileparstmp[256] = ""; FILE *fppars = NULL, *fpparstmp = NULL; UINT4 nonfixed = 1; outfile = ppt1->value; /* open file for printing out list of parameter names - this should already exist */ sprintf(outfilepars, "%s_params.txt", outfile); if( (fppars = fopen(outfilepars, "r")) == NULL ){ XLAL_ERROR_VOID(XLAL_EIO, "Error... cannot open parameter name output file %s.", outfilepars); } /* read in the parameter names and remove the "model" value */ sprintf(outfileparstmp, "%s_params.txt_tmp", outfile); if( (fpparstmp = fopen(outfileparstmp, "w")) == NULL ){ XLAL_ERROR_VOID(XLAL_EIO, "Error... cannot open parameter name output file %s.", outfileparstmp); } if ( LALInferenceGetProcParamVal( runState->commandLine, "--output-all-params" ) ){ nonfixed = 0; } CHAR v[128] = ""; while( fscanf(fppars, "%s", v) != EOF ){ /* if outputing only non-fixed values then only re-output names of those non-fixed things */ if ( nonfixed ){ if ( LALInferenceCheckVariable( threadState->currentParams, v ) ){ if ( LALInferenceGetVariableVaryType( threadState->currentParams, v ) != LALINFERENCE_PARAM_FIXED ){ fprintf(fpparstmp, "%s\t", v); } } } else{ /* re-output everything to a temporary file */ fprintf(fpparstmp, "%s\t", v); } } fclose(fppars); fclose(fpparstmp); /* move the temporary file name to the standard outfile_param name */ rename( outfileparstmp, outfilepars ); /* gzip the output file if required */ if( LALInferenceGetProcParamVal( runState->commandLine, "--gzip" ) ){ if ( XLALGzipTextFile( outfile ) != XLAL_SUCCESS ){ XLAL_PRINT_ERROR( "Error... Could not gzip the output file!\n" ); XLAL_ERROR_VOID( XLAL_EIO ); } } } /* if we have XML enabled */ #ifdef HAVE_LIBLALXML ProcessParamsTable *ppt2 = LALInferenceGetProcParamVal( runState->commandLine, "--outxml" ); if(!ppt2){ ppt2=LALInferenceGetProcParamVal(runState->commandLine,"--outXML"); } CHAR *outVOTable = NULL; if ( !ppt2 && !ppt1 ){ XLAL_PRINT_ERROR("Must specify either --outfile or --outXML\n"); XLAL_ERROR_VOID( XLAL_EIO ); } if( ppt2 ){ outVOTable = ppt2->value; /* gzip the output file if required */ if( LALInferenceGetProcParamVal( runState->commandLine, "--gzip" ) ){ if ( XLALGzipTextFile( outVOTable ) != XLAL_SUCCESS ){ XLAL_ERROR_VOID( XLAL_EIO, "Error... Could not gzip the output file!" ); } } } #else if ( !ppt1 ){ XLAL_ERROR_VOID( XLAL_EIO, "Error... --outfile not defined!" ); } #endif return; }
/** * @brief Advance a LALFrStream stream to the beginning of the next frame * @details * The position of a LALFrStream is advanced so that the next read will * be at the next frame. If the stream is at the end, the #LAL_FR_STREAM_END * bit of the LALFrStreamState state is set, and the routine returns the * return code 1. If there is a gap in the data before the next frame, * the #LAL_FR_STREAM_GAP bit of the LALFrStreamState state is set, and the * routine returns the return code 2. If, however, the * #LAL_FR_STREAM_IGNOREGAP_MODE bit is not set in the LALFrStreamMode mode * then the routine produces an error if a gap is encountered. * @param stream Pointer to a #LALFrStream structure. * @retval 2 Gap in the data is encountered. * @retval 1 End of stream encountered. * @retval 0 Normal success. * @retval <0 Failure. */ int XLALFrStreamNext(LALFrStream * stream) { /* timing accuracy: tenth of a sample interval for a 16kHz fast channel */ const INT8 tacc = (INT8) floor(0.1 * 1e9 / 16384.0); const char *url1; const char *url2; int pos1; int pos2; INT8 texp = 0; INT8 tact; if (stream->state & LAL_FR_STREAM_END) return 1; /* end code */ /* turn off gap bit */ stream->state &= ~LAL_FR_STREAM_GAP; url2 = url1 = stream->cache->list[stream->fnum].url; pos2 = pos1 = stream->pos; /* open a new file if necessary */ if (!stream->file) { if (stream->fnum >= stream->cache->length) { stream->state |= LAL_FR_STREAM_END; return 1; } if (XLALFrStreamFileOpen(stream, stream->fnum) < 0) XLAL_ERROR(XLAL_EFUNC); } if (stream->file) { INT4 nFrame = XLALFrFileQueryNFrame(stream->file); if (stream->pos < nFrame) { LIGOTimeGPS gpstime; XLALGPSToINT8NS(XLALFrFileQueryGTime(&gpstime, stream->file, stream->pos)); texp = XLALGPSToINT8NS(XLALGPSAdd(&gpstime, XLALFrFileQueryDt(stream->file, stream->pos))); ++stream->pos; } if (stream->pos >= nFrame) { XLALFrStreamFileClose(stream); ++stream->fnum; } pos2 = stream->pos; } /* open a new file if necessary */ if (!stream->file) { if (stream->fnum >= stream->cache->length) { stream->state |= LAL_FR_STREAM_END; return 1; } if (XLALFrStreamFileOpen(stream, stream->fnum) < 0) XLAL_ERROR(XLAL_EFUNC); url2 = stream->cache->list[stream->fnum].url; pos2 = stream->pos; } /* compute actual start time of this new frame */ tact = XLALGPSToINT8NS(XLALFrFileQueryGTime(&stream->epoch, stream->file, stream->pos)); /* INT8 is platform dependent, cast to long long for llabs() call */ if (llabs((long long)(texp - tact)) > tacc) { /* there is a gap */ stream->state |= LAL_FR_STREAM_GAP; if (stream->mode & LAL_FR_STREAM_GAPINFO_MODE) { XLAL_PRINT_INFO("Gap in frame data between times %.6f and %.6f", 1e-9 * texp, 1e-9 * tact); } if (!(stream->mode & LAL_FR_STREAM_IGNOREGAP_MODE)) { XLAL_PRINT_ERROR("Gap in frame data"); XLAL_PRINT_ERROR("Time %.6f is end of frame %d of file %s", 1e-9 * texp, pos1, url1); XLAL_PRINT_ERROR("Time %.6f is start of frame %d of file %s", 1e-9 * tact, pos2, url2); XLAL_ERROR(XLAL_ETIME); } return 2; /* gap code */ } return 0; }