int LALInferenceInspiralPriorTest(void)
{
	TEST_HEADER();

	int errnum;

	REAL8 result;
	LALInferenceRunState *runState = XLALCalloc(1, sizeof(LALInferenceRunState));
    LALInferenceThreadState *thread = XLALCalloc(1, sizeof(LALInferenceThreadState));
    runState->threads = XLALCalloc(1, sizeof(LALInferenceThreadState*));
    runState->threads[0] = thread;
	LALInferenceVariables *params = XLALCalloc(1, sizeof(LALInferenceVariables));
	LALInferenceVariables *priorArgs = XLALCalloc(1, sizeof(LALInferenceVariables));

	// Standard null reference check.
	int failed = 1;
	runState->priorArgs = NULL;
    thread->model = NULL;
	XLAL_TRY(result = LALInferenceInspiralPrior(runState, params, thread->model), errnum);
	failed &= !XLAL_IS_REAL8_FAIL_NAN(result) || errnum != XLAL_EFAULT;
	runState->priorArgs = priorArgs;
	XLAL_TRY(result = LALInferenceInspiralPrior(NULL, params, thread->model), errnum);
	failed &= !XLAL_IS_REAL8_FAIL_NAN(result) || errnum != XLAL_EFAULT;
	XLAL_TRY(result = LALInferenceInspiralPrior(runState, NULL, thread->model), errnum);
	failed &= !XLAL_IS_REAL8_FAIL_NAN(result) || errnum != XLAL_EFAULT;
	if (failed)
		TEST_FAIL("Null reference check failed.");

	// Set up parameters.
	REAL8 value = 0;
	LALInferenceAddVariable(params, "logdistance", &value, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_LINEAR);
	LALInferenceAddVariable(params, "distance", &value, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_LINEAR);
	LALInferenceAddVariable(params, "inclination", &value, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_CIRCULAR);
	LALInferenceAddVariable(params, "rightascension", &value, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_CIRCULAR);
	LALInferenceAddVariable(params, "declination", &value, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_CIRCULAR);
	LALInferenceAddVariable(params, "theta_spin1", &value, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_CIRCULAR);
	LALInferenceAddVariable(params, "theta_spin2", &value, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_CIRCULAR);
	LALInferenceAddVariable(params, "logmc", &value, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_LINEAR);
	LALInferenceAddVariable(params, "chirpmass", &value, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_LINEAR);
	/*LALInferenceAddVariable(params, "q", &value, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_LINEAR);*/
	LALInferenceAddVariable(params, "eta", &value, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_LINEAR);
	
	REAL8 min, max;
	min = 1.0; max = 100.0;
	LALInferenceAddMinMaxPrior(priorArgs, "distance", &min, &max, LALINFERENCE_REAL8_t);
	min = log(min); max = log(max);
	LALInferenceAddMinMaxPrior(priorArgs, "logdistance", &min, &max, LALINFERENCE_REAL8_t);
	min = 1.0; max = 20.5;
	LALInferenceAddMinMaxPrior(priorArgs, "chirpmass", &min, &max, LALINFERENCE_REAL8_t);
	min = log(min); max = log(max);
	LALInferenceAddMinMaxPrior(priorArgs, "logmc", &min, &max, LALINFERENCE_REAL8_t);
	min = 1.0; max = 30.0;
	LALInferenceAddMinMaxPrior(priorArgs, "mass1", &min, &max, LALINFERENCE_REAL8_t);
	LALInferenceAddMinMaxPrior(priorArgs, "mass2", &min, &max, LALINFERENCE_REAL8_t);
	/*max *= 2;*/
	/*LALInferenceAddVariable(priorArgs, "MTotMax", &max, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_FIXED);*/
	min = -LAL_PI; max = LAL_PI;
	LALInferenceAddMinMaxPrior(priorArgs, "inclination", &min, &max, LALINFERENCE_REAL8_t);
	min = 0; max = LAL_TWOPI;
	LALInferenceAddMinMaxPrior(priorArgs, "rightascension", &min, &max, LALINFERENCE_REAL8_t);
	min = -LAL_PI / 2.0; max = LAL_PI / 2.0;
	LALInferenceAddMinMaxPrior(priorArgs, "declination", &min, &max, LALINFERENCE_REAL8_t);
	min = -LAL_PI / 2.0; max = LAL_PI / 2.0;
	LALInferenceAddMinMaxPrior(priorArgs, "theta_spin1", &min, &max, LALINFERENCE_REAL8_t);
	min = -LAL_PI / 2.0; max = LAL_PI / 2.0;
	LALInferenceAddMinMaxPrior(priorArgs, "theta_spin2", &min, &max, LALINFERENCE_REAL8_t);
	min = 0.01; max = 0.25;
	LALInferenceAddMinMaxPrior(priorArgs, "eta", &min, &max, LALINFERENCE_REAL8_t);

	// Pick a random point in the non-zero region of the parameter space.
	gsl_rng *rng = gsl_rng_alloc(gsl_rng_default);
	gsl_rng_set(rng, 0);
	LALInferenceDrawFromPrior(params, priorArgs, rng);

	// Check that we get a finite log prior.
	XLAL_TRY(result = LALInferenceInspiralPrior(runState, params, thread->model), errnum);

	if (XLAL_IS_REAL8_FAIL_NAN(result) || errnum != XLAL_SUCCESS)
	{
		TEST_FAIL("Could not generate inspiral prior; XLAL error: %s", XLALErrorString(errnum));
	}
	else if (result == -DBL_MAX)
	{
		TEST_FAIL("Parameter configuration within specified min/max bounds for each parameter gave zero prior.");
	}

	// Now set a parameter outside its bounds and see what happens.
	LALInferenceGetMinMaxPrior(priorArgs, "distance", &min, &max);
	value = max + (max - min) / 2;
	LALInferenceSetVariable(params, "distance", &value);
	XLAL_TRY(result = LALInferenceInspiralPrior(runState, params, thread->model), errnum);
	if (XLAL_IS_REAL8_FAIL_NAN(result) || errnum != XLAL_SUCCESS)
	{
		TEST_FAIL("Could not generate inspiral prior; XLAL error: %s", XLALErrorString(errnum));
	}
	else if (result != -DBL_MAX)
	{
		TEST_FAIL("Distance %f is outside [%f,%f] but prior is non-zero.", value, min, max);
	}

	// Try another configuration; this time set m1 and m2 such that one is *outside* its bounds,
	// but the chirp mass and symmetric mass ratio are still OK; this should be picked up and a
	// zero prior returned.
	LALInferenceDrawFromPrior(params, priorArgs, rng);
	LALInferenceGetMinMaxPrior(priorArgs, "mass1", &min, &max);
	REAL8 m2 = 0.5;
	REAL8 m1 = 3.82;
	REAL8 eta = m1 * m2 / pow(m1 + m2, 2);
	LALInferenceSetVariable(params, "eta", &eta);
	REAL8 Mc = pow(m1 * m2, 3.0 / 5.0) / pow(m1 + m2, 1.0 / 5.0);
	LALInferenceSetVariable(params, "chirpmass", &Mc);
	REAL8 logMc = log(Mc);
	LALInferenceSetVariable(params, "logmc", &logMc);
	XLAL_TRY(result = LALInferenceInspiralPrior(runState, params, thread->model), errnum);
	if (XLAL_IS_REAL8_FAIL_NAN(result) || errnum != XLAL_SUCCESS)
	{
		TEST_FAIL("Could not generate inspiral prior; XLAL error: %s", XLALErrorString(errnum));
	}
	else if (result != -DBL_MAX)
	{
		TEST_FAIL("Mass ratio %f and chirp mass %f define masses outside bounds [%f,%f], but prior is non-zero.", eta, Mc, min, max);
	}

	TEST_FOOTER();
}
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;
}
int LALInferenceRotateInitialPhaseTest(void)
{
	TEST_HEADER();

	int errnum;

	// A basic null reference check.
	XLAL_TRY(LALInferenceRotateInitialPhase(NULL), errnum);
	if (errnum != XLAL_EFAULT)
		TEST_FAIL("Null reference check failed.");

	// Construct a variable list containing phi0 and psi.
	REAL8 psi = 0;
	REAL8 phi0 = 0;
	REAL8 phi0_2 = 0;
	LALInferenceVariables *variables = XLALCalloc(1, sizeof(LALInferenceVariables));
	LALInferenceAddVariable(variables, "psi", &psi, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_CIRCULAR);
	LALInferenceAddVariable(variables, "phi0", &phi0, LALINFERENCE_REAL8_t, LALINFERENCE_PARAM_CIRCULAR);

	// Check that if psi is in [0,2pi], phi0 remains unchanged.
	psi = LAL_PI;
	phi0 = 0;
	LALInferenceSetVariable(variables, "psi", &psi);
	LALInferenceSetVariable(variables, "phi0", &phi0);
	LALInferenceRotateInitialPhase(variables);
	phi0_2 = *(REAL8 *)LALInferenceGetVariable(variables, "phi0");
	if (!compareFloats(phi0_2, phi0, EPSILON))
		TEST_FAIL("Psi in [0,2pi] but phi0 has changed!");

	// Check that if psi is outside [0,2pi], phi0 is rotated.
	psi = -2 * LAL_TWOPI;
	phi0 = 0;
	LALInferenceSetVariable(variables, "psi", &psi);
	LALInferenceSetVariable(variables, "phi0", &phi0);
	LALInferenceRotateInitialPhase(variables);
	phi0_2 = *(REAL8 *)LALInferenceGetVariable(variables, "phi0");
	if (!compareFloats(phi0_2, phi0 + LAL_PI, EPSILON))
		TEST_FAIL("Psi outside [0,2pi] but phi0 not rotated by 2pi!");
	psi = 2 * LAL_TWOPI;
	phi0 = 0;
	LALInferenceSetVariable(variables, "psi", &psi);
	LALInferenceSetVariable(variables, "phi0", &phi0);
	LALInferenceRotateInitialPhase(variables);
	phi0_2 = *(REAL8 *)LALInferenceGetVariable(variables, "phi0");
	if (!compareFloats(phi0_2, phi0 + LAL_PI, EPSILON))
		TEST_FAIL("Psi outside [0,2pi] but phi0 not rotated by 2pi!");

	// Check boundary cases: if psi=0 or psi=2pi, phi0 shouldn't be rotated.
	psi = 0;
	phi0 = 0;
	LALInferenceSetVariable(variables, "psi", &psi);
	LALInferenceSetVariable(variables, "phi0", &phi0);
	LALInferenceRotateInitialPhase(variables);
	phi0_2 = *(REAL8 *)LALInferenceGetVariable(variables, "phi0");
	if (phi0 != phi0_2) // Should be exact.)
		TEST_FAIL("Psi on boundary of [0,2pi] but phi0 has changed!");
	psi = LAL_TWOPI;
	phi0 = 0;
	LALInferenceSetVariable(variables, "psi", &psi);
	LALInferenceSetVariable(variables, "phi0", &phi0);
	LALInferenceRotateInitialPhase(variables);
	phi0_2 = *(REAL8 *)LALInferenceGetVariable(variables, "phi0");
	if (phi0 != phi0_2) // Should be exact.)
		TEST_FAIL("Psi on boundary of [0,2pi] but phi0 has changed!");

	XLALFree(variables);
	TEST_FOOTER();
}