/**
 * Function to calculate the second derivative of the Hamiltonian.
* The derivatives are taken with respect to indices idx1, idx2    */
static REAL8
XLALCalculateSphHamiltonianDeriv2Prec(
				  const int idx1,	/**<< Derivative w.r.t. index 1 */
				  const int idx2,	/**<< Derivative w.r.t. index 2 */
				  const REAL8 values[],	/**<< Dynamical variables in spherical coordinates */
				  SpinEOBParams * params	/**<< Spin EOB Parameters */
)
{

	static const REAL8 STEP_SIZE = 1.0e-5;

	REAL8		result;
	REAL8 UNUSED	absErr;

	HcapSphDeriv2Params dParams;

	gsl_function	F;
	INT4 UNUSED	gslStatus;

	dParams.sphValues = values;
	dParams.varyParam1 = idx1;
	dParams.varyParam2 = idx2;
	dParams.params = params;

	/*
	 * XLAL_PRINT_INFO( " In second deriv function: values\n" ); for ( int i = 0;
	 * i < 12; i++ ) { XLAL_PRINT_INFO( "%.16e ", values[i] ); } XLAL_PRINT_INFO( "\n" );
	 */
	F.function = GSLSpinHamiltonianDerivWrapperPrec;
	F.params = &dParams;

	/* GSL seemed to give weird answers - try my own code */
	/*
	 * result = GSLSpinHamiltonianDerivWrapperPrec( values[idx1] + STEP_SIZE,
	 * &dParams ) - GSLSpinHamiltonianDerivWrapperPrec( values[idx1] -
	 * STEP_SIZE, &dParams ); XLAL_PRINT_INFO( "%.16e - %.16e / 2h\n",
	 * GSLSpinHamiltonianDerivWrapperPrec( values[idx1] + STEP_SIZE, &dParams
	 * ), GSLSpinHamiltonianDerivWrapperPrec( values[idx1] - STEP_SIZE,
	 * &dParams ) );
	 *
	 * result = result / ( 2.*STEP_SIZE );
	 */

	XLAL_CALLGSL(gslStatus = gsl_deriv_central(&F, values[idx1],
					      STEP_SIZE, &result, &absErr));

	if (gslStatus != GSL_SUCCESS) {
		XLALPrintError("XLAL Error %s - Failure in GSL function\n", __func__);
		XLAL_ERROR_REAL8(XLAL_EDOM);
	}
	//XLAL_PRINT_INFO("Second deriv abs err = %.16e\n", absErr);

	//XLAL_PRINT_INFO("RESULT = %.16e\n", result);
	return result;
}
/*
-------------------------------------------------------------------------
 * This function frees the memory allocated using XLALInitFContactWorkSpace.
 * ------------------------------------------------------------------------*/
void XLALFreeFContactWorkSpace( fContactWorkSpace *workSpace )
{
  if ( !workSpace )
    XLAL_ERROR_VOID( XLAL_EFAULT );

  /* Free all the allocated memory */
  if (workSpace->tmpA) XLAL_CALLGSL( gsl_matrix_free( workSpace->tmpA ));
  if (workSpace->tmpB) XLAL_CALLGSL( gsl_matrix_free( workSpace->tmpB ));
  if (workSpace->C)    XLAL_CALLGSL( gsl_matrix_free( workSpace->C ));
  if (workSpace->p1)   XLAL_CALLGSL( gsl_permutation_free( workSpace->p1 ));
  if (workSpace->tmpV) XLAL_CALLGSL( gsl_vector_free( workSpace->tmpV ));
  if (workSpace->r_AB) XLAL_CALLGSL( gsl_vector_free( workSpace->r_AB ));
  if (workSpace->s)    XLAL_CALLGSL( gsl_min_fminimizer_free( workSpace->s ));

  LALFree( workSpace );
  return;
}
/*------------------------------------------------------------------------------------------
 *
 *    Definitions of functions (only one in this file, so no prototypes needed).
 *
 *------------------------------------------------------------------------------------------
 */
static int SEOBNRv3OptimizedInterpolatorGeneral(
                REAL8 * yin, /**<< Data to be interpolated; time first */
                REAL8 tinit, /**<< time at which to begin interpolating */
                REAL8 deltat, /**<< Spacing between interpolated times */
                UINT4 num_input_times, /**<< The number of input times */
                REAL8Array ** yout, /**<< Interpolation output */
                size_t dim /**<< Number of quantities interpolated (e.g. if yin = {t,x,y,z} then dim 3) */
                )
{
    int errnum = 0;

    /* needed for the final interpolation */
    gsl_spline *interp = NULL;
    gsl_interp_accel *accel = NULL;
    int outputlen = 0;
    REAL8Array *output = NULL;
    REAL8 *times, *vector;      /* aliases */

    /* note: for speed, this replaces the single CALLGSL wrapper applied before each GSL call */
    interp = gsl_spline_alloc(gsl_interp_cspline, num_input_times);
    accel = gsl_interp_accel_alloc();

    outputlen = (int)(yin[num_input_times-1] / deltat) + 1;

    output = XLALCreateREAL8ArrayL(2, dim+1, outputlen);/* Original (dim+1)*/

    if (!interp || !accel || !output) {
      errnum = XLAL_ENOMEM;   /* ouch again, ran out of memory */
      if (output)
        XLALDestroyREAL8Array(output);
      outputlen = 0;
      goto bail_out;
    }

    /* make an array of times */
    times = output->data;
    for (int j = 0; j < outputlen; j++)
      times[j] = tinit + deltat * j;

    /* interpolate! */
    for (unsigned int i = 1; i <= dim; i++) { /* Original (dim)  */
     //gsl_spline_init(interp, &yin->data[0], &yin->data[num_input_times * i], num_input_times + 1);
     gsl_spline_init(interp, yin, &(yin[num_input_times * i]), num_input_times);

      vector = output->data + outputlen * i;
      unsigned int index_old=0;
      double x_lo_old=0,y_lo_old=0,b_i_old=0,c_i_old=0,d_i_old=0;
      for (int j = 0; j < outputlen; j++) {
        optimized_gsl_spline_eval_e(interp,times[j],accel, &(vector[j]),&index_old,&x_lo_old,&y_lo_old,&b_i_old,&c_i_old,&d_i_old);
      }
    }

    /* deallocate stuff and return */
  bail_out:

    if (interp)
        XLAL_CALLGSL(gsl_spline_free(interp));
    if (accel)
        XLAL_CALLGSL(gsl_interp_accel_free(accel));

    if (errnum)
        XLAL_ERROR(errnum);

    *yout = output;
    return outputlen;
}
/**
 * This function is largely based on/copied from
 * XLALAdaptiveRungeKutta4(), which exists inside the
 * lal/src/utilities/LALAdaptiveRungeKutta4.c file
 * subroutine. It reads in an array of timeseries that
 * contain data *not* evenly spaced in time
 * and performs cubic spline interpolations to resample
 * the data to uniform time sampling. Interpolations use
 * GSL's built-in routines, which recompute interpolation
 * coefficients each time the itnerpolator is called.
 * This can be extremely inefficient; in case of SEOBNRv4,
 * first data points exist at very large dt, and
 * interp. coefficients might be needlessly recomputed
 * 10,000+ times or more. We also made the optimization
 * that assumes the data are sampled at points monotone
 * in time.
 * tinit and deltat specify the desired initial time and
 *   time spacing for the output interpolated data,
 * num_input_times denotes the number of points yin arrays
 *   are sampled in time.
 * yout is the output array.
 */
UNUSED static int
SEOBNRv2OptimizedInterpolatorNoAmpPhase (REAL8Array * yin, REAL8 tinit,
					 REAL8 deltat, UINT4 num_input_times,
					 REAL8Array ** yout)
{
  int errnum = 0;

  /* needed for the final interpolation */
  gsl_spline *interp = NULL;
  gsl_interp_accel *accel = NULL;
  int outputlen = 0;
  REAL8Array *output = NULL;
  REAL8 *times, *vector;	/* aliases */

  /* needed for the integration */
  size_t dim = 4;
  
  interp = gsl_spline_alloc (gsl_interp_cspline, num_input_times);
  accel = gsl_interp_accel_alloc ();

  outputlen = (int) (yin->data[num_input_times - 1] / deltat) + 1;
  output = XLALCreateREAL8ArrayL (2, dim + 1, outputlen);	/* Only dim + 1 rather than dim + 3 since we're not adding amp & phase */

  if (!interp || !accel || !output)
    {
      errnum = XLAL_ENOMEM;	/* ouch again, ran out of memory */
      if (output)
	XLALDestroyREAL8Array (output);
      outputlen = 0;
      goto bail_out;
    }

  /* make an array of times */
  times = output->data;
  for (int j = 0; j < outputlen; j++)
    times[j] = tinit + deltat * j;

  /* interpolate! */
  for (unsigned int i = 1; i <= dim; i++)
    {				/* only up to dim (4) because we are not interpolating amplitude and phase */
      //gsl_spline_init(interp, &yin->data[0], &yin->data[num_input_times * i], num_input_times + 1);
      gsl_spline_init (interp, &yin->data[0], &yin->data[num_input_times * i],
		       num_input_times);

      vector = output->data + outputlen * i;
      unsigned int index_old = 0;
      double x_lo_old = 0, y_lo_old = 0, b_i_old = 0, c_i_old = 0, d_i_old =
	0;
      for (int j = 0; j < outputlen; j++)
	{
	  optimized_gsl_spline_eval_e (interp, times[j], accel, &(vector[j]),
				       &index_old, &x_lo_old, &y_lo_old,
				       &b_i_old, &c_i_old, &d_i_old);
	}
    }

  /* deallocate stuff and return */
bail_out:

  if (interp)
    XLAL_CALLGSL (gsl_spline_free (interp));
  if (accel)
    XLAL_CALLGSL (gsl_interp_accel_free (accel));

  if (errnum)
    XLAL_ERROR (errnum);

  *yout = output;
  return outputlen;
}
int main(void) {
  gsl_matrix *TS; /* the training set of real waveforms */
  gsl_matrix_complex *cTS; /* the training set of complex waveforms */

  size_t TSsize;  /* the size of the training set (number of waveforms) */
  size_t wl;      /* the length of each waveform */
  size_t k = 0, j = 0, i = 0, nbases = 0, cnbases = 0;

  REAL8 *RB = NULL;        /* the real reduced basis set */
  COMPLEX16 *cRB = NULL;   /* the complex reduced basis set */
  LALInferenceREALROQInterpolant *interp = NULL;
  LALInferenceCOMPLEXROQInterpolant *cinterp = NULL;

  gsl_vector *freqs;

  double tolerance = TOLERANCE; /* tolerance for reduced basis generation loop */

  TSsize = TSSIZE;
  wl = WL;

  /* allocate memory for training set */
  TS = gsl_matrix_calloc(TSsize, wl);
  cTS = gsl_matrix_complex_calloc(TSsize, wl);

  /* the waveform model is just a simple chirp so set up chirp mass range for training set */
  double fmin0 = 48, fmax0 = 256, f0 = 0., m0 = 0.;
  double df = (fmax0-fmin0)/(wl-1.); /* model time steps */
  freqs = gsl_vector_alloc(wl); /* times at which to calculate the model */
  double Mcmax = 2., Mcmin = 1.5, Mc = 0.;

  gsl_vector_view fweights = gsl_vector_view_array(&df, 1);

  /* set up training sets (one real and one complex) */
  for ( k=0; k < TSsize; k++ ){
    Mc = pow(pow(Mcmin, 5./3.) + (double)k*(pow(Mcmax, 5./3.)-pow(Mcmin, 5./3.))/((double)TSsize-1), 3./5.);

    for ( j=0; j < wl; j++ ){
      f0 = fmin0 + (double)j*(fmax0-fmin0)/((double)wl-1.);

      gsl_complex gctmp;
      COMPLEX16 ctmp;
      m0 = real_model(f0, Mc);
      ctmp = imag_model(f0, Mc);
      GSL_SET_COMPLEX(&gctmp, creal(ctmp), cimag(ctmp));
      gsl_vector_set(freqs, j, f0);
      gsl_matrix_set(TS, k, j, m0);
      gsl_matrix_complex_set(cTS, k, j, gctmp);
    }
  }

  /* create reduced orthonormal basis from training set */
  if ( (RB = LALInferenceGenerateREAL8OrthonormalBasis(&fweights.vector, tolerance, TS, &nbases)) == NULL){
    fprintf(stderr, "Error... problem producing basis\n");
    return 1;
  }

  if ( (cRB = LALInferenceGenerateCOMPLEX16OrthonormalBasis(&fweights.vector, tolerance, cTS, &cnbases)) == NULL){
    fprintf(stderr, "Error... problem producing basis\n");
    return 1;
  }

  /* free the training set */
  gsl_matrix_free(TS);
  gsl_matrix_complex_free(cTS);

  gsl_matrix_view RBview = gsl_matrix_view_array(RB, nbases, wl);
  gsl_matrix_complex_view cRBview = gsl_matrix_complex_view_array((double*)cRB, cnbases, wl);

  fprintf(stderr, "No. nodes (real)  = %zu, %zu x %zu\n", nbases, RBview.matrix.size1, RBview.matrix.size2);
  fprintf(stderr, "No. nodes (complex)  = %zu, %zu x %zu\n", cnbases, cRBview.matrix.size1, cRBview.matrix.size2);

  /* get the interpolant */
  interp = LALInferenceGenerateREALROQInterpolant(&RBview.matrix);
  cinterp = LALInferenceGenerateCOMPLEXROQInterpolant(&cRBview.matrix);

  /* free the reduced basis */
  XLALFree(RB);
  XLALFree(cRB);

  /* now get the terms for the likelihood with and without the reduced order quadrature
   * and do some timing tests */

  /* create the model dot model weights */
  REAL8 varval = 1.;
  gsl_vector_view vars = gsl_vector_view_array(&varval, 1);

  gsl_matrix *mmw = LALInferenceGenerateREALModelModelWeights(interp->B, &vars.vector);
  gsl_matrix_complex *cmmw = LALInferenceGenerateCOMPLEXModelModelWeights(cinterp->B, &vars.vector);

  /* let's create some Gaussian random data */
  const gsl_rng_type *T;
  gsl_rng *r;

  gsl_rng_env_setup();

  T = gsl_rng_default;
  r = gsl_rng_alloc(T);

  REAL8 *data = XLALCalloc(wl, sizeof(REAL8));
  COMPLEX16 *cdata = XLALCalloc(wl, sizeof(COMPLEX16));
  for ( i=0; i<wl; i++ ){
    data[i] = gsl_ran_gaussian(r, 1.0);                               /* real data */
    cdata[i] = gsl_ran_gaussian(r, 1.0) + I*gsl_ran_gaussian(r, 1.0); /* complex data */
  }

  /* create the data dot model weights */
  gsl_vector_view dataview = gsl_vector_view_array(data, wl);
  gsl_vector *dmw = LALInferenceGenerateREAL8DataModelWeights(interp->B, &dataview.vector, &vars.vector);

  gsl_vector_complex_view cdataview = gsl_vector_complex_view_array((double*)cdata, wl);
  gsl_vector_complex *cdmw = LALInferenceGenerateCOMPLEX16DataModelWeights(cinterp->B, &cdataview.vector, &vars.vector);

  /* pick a chirp mass and generate a model to compare likelihoods */
  double randMc = 1.873; /* a random frequency to create a model */

  gsl_vector *modelfull = gsl_vector_alloc(wl);
  gsl_vector *modelreduced = gsl_vector_alloc(nbases);
  gsl_vector_complex *cmodelfull = gsl_vector_complex_alloc(wl);
  gsl_vector_complex *cmodelreduced = gsl_vector_complex_alloc(cnbases);

  /* create models */
  for ( i=0; i<wl; i++ ){
    /* models at all frequencies */
    gsl_vector_set(modelfull, i, real_model(gsl_vector_get(freqs, i), randMc));

    COMPLEX16 cval = imag_model(gsl_vector_get(freqs, i), randMc);
    gsl_complex gcval;
    GSL_SET_COMPLEX(&gcval, creal(cval), cimag(cval));
    gsl_vector_complex_set(cmodelfull, i, gcval);
  }

  /* models at interpolant nodes */
  for ( i=0; i<nbases; i++ ){ /* real model */
    gsl_vector_set(modelreduced, i, real_model(gsl_vector_get(freqs, interp->nodes[i]), randMc));
  }
  for ( i=0; i<cnbases; i++ ){ /* complex model */
    COMPLEX16 cval = imag_model(gsl_vector_get(freqs, cinterp->nodes[i]), randMc);
    gsl_complex gcval;
    GSL_SET_COMPLEX(&gcval, creal(cval), cimag(cval));
    gsl_vector_complex_set(cmodelreduced, i, gcval);
  }

  /* timing variables */
  struct timeval t1, t2, t3, t4;
  double dt1, dt2;

  /* start with the real model */
  /* get the model model term with the full model */
  REAL8 mmfull, mmred;
  gettimeofday(&t1, NULL);
  XLAL_CALLGSL( gsl_blas_ddot(modelfull, modelfull, &mmfull) );        /* real model */
  gettimeofday(&t2, NULL);

  /* now get it with the reduced order quadrature */
  gettimeofday(&t3, NULL);
  mmred = LALInferenceROQREAL8ModelDotModel(mmw, modelreduced);
  gettimeofday(&t4, NULL);

  dt1 = (double)((t2.tv_sec + t2.tv_usec*1.e-6) - (t1.tv_sec + t1.tv_usec*1.e-6));
  dt2 = (double)((t4.tv_sec + t4.tv_usec*1.e-6) - (t3.tv_sec + t3.tv_usec*1.e-6));
  fprintf(stderr, "Real Signal:\n - M dot M (full) = %le [%.9lf s], M dot M (reduced) = %le [%.9lf s], time ratio = %lf\n", mmfull, dt1, mmred, dt2, dt1/dt2);

  /* get the data model term with the full model */
  REAL8 dmfull, dmred;
  gettimeofday(&t1, NULL);
  XLAL_CALLGSL( gsl_blas_ddot(&dataview.vector, modelfull, &dmfull) );
  gettimeofday(&t2, NULL);

  /* now get it with the reduced order quadrature */
  gettimeofday(&t3, NULL);
  dmred = LALInferenceROQREAL8DataDotModel(dmw, modelreduced);
  gettimeofday(&t4, NULL);

  dt1 = (double)((t2.tv_sec + t2.tv_usec*1.e-6) - (t1.tv_sec + t1.tv_usec*1.e-6));
  dt2 = (double)((t4.tv_sec + t4.tv_usec*1.e-6) - (t3.tv_sec + t3.tv_usec*1.e-6));
  fprintf(stderr, " - D dot M (full) = %le [%.9lf s], D dot M (reduced) = %le [%.9lf s], time ratio = %lf\n", dmfull, dt1, dmred, dt2, dt1/dt2);

  /* check difference in log likelihoods */
  double Lfull, Lred, Lfrac;

  Lfull = mmfull - 2.*dmfull;
  Lred = mmred - 2.*dmred;
  Lfrac = 100.*fabs(Lfull-Lred)/fabs(Lfull); /* fractional log likelihood difference (in %) */

  fprintf(stderr, " - Fractional difference in log likelihoods = %lf%%\n", Lfrac);

  XLALFree(data);
  gsl_vector_free(modelfull);
  gsl_vector_free(modelreduced);
  gsl_matrix_free(mmw);
  gsl_vector_free(dmw);

  /* check log likelihood difference is within tolerance */
  if ( Lfrac > LTOL ) { return 1; }

  /* now do the same with the complex model */
  /* get the model model term with the full model */
  COMPLEX16 cmmred, cmmfull;
  gsl_complex cmmfulltmp;
  gettimeofday(&t1, NULL);
  XLAL_CALLGSL( gsl_blas_zdotc(cmodelfull, cmodelfull, &cmmfulltmp) ); /* complex model */
  cmmfull = GSL_REAL(cmmfulltmp) + I*GSL_IMAG(cmmfulltmp);
  gettimeofday(&t2, NULL);

  gettimeofday(&t3, NULL);
  cmmred = LALInferenceROQCOMPLEX16ModelDotModel(cmmw, cmodelreduced);
  gettimeofday(&t4, NULL);

  dt1 = (double)((t2.tv_sec + t2.tv_usec*1.e-6) - (t1.tv_sec + t1.tv_usec*1.e-6));
  dt2 = (double)((t4.tv_sec + t4.tv_usec*1.e-6) - (t3.tv_sec + t3.tv_usec*1.e-6));
  fprintf(stderr, "Complex Signal:\n - M dot M (full) = %le [%.9lf s], M dot M (reduced) = %le [%.9lf s], time ratio = %lf\n", creal(cmmfull), dt1, creal(cmmred), dt2, dt1/dt2);

  COMPLEX16 cdmfull, cdmred;
  gsl_complex cdmfulltmp;
  gettimeofday(&t1, NULL);
  XLAL_CALLGSL( gsl_blas_zdotc(&cdataview.vector, cmodelfull, &cdmfulltmp) );
  cdmfull = GSL_REAL(cdmfulltmp) + I*GSL_IMAG(cdmfulltmp);
  gettimeofday(&t2, NULL);

  gettimeofday(&t3, NULL);
  cdmred = LALInferenceROQCOMPLEX16DataDotModel(cdmw, cmodelreduced);
  gettimeofday(&t4, NULL);

  dt1 = (double)((t2.tv_sec + t2.tv_usec*1.e-6) - (t1.tv_sec + t1.tv_usec*1.e-6));
  dt2 = (double)((t4.tv_sec + t4.tv_usec*1.e-6) - (t3.tv_sec + t3.tv_usec*1.e-6));
  fprintf(stderr, " - D dot M (full) = %le [%.9lf s], D dot M (reduced) = %le [%.9lf s], time ratio = %lf\n", creal(cdmfull), dt1, creal(cdmred), dt2, dt1/dt2);

  /* check difference in log likelihoods */
  Lfull = creal(cmmfull) - 2.*creal(cdmfull);
  Lred = creal(cmmred) - 2.*creal(cdmred);
  Lfrac = 100.*fabs(Lfull-Lred)/fabs(Lfull); /* fractional log likelihood difference (in %) */

  fprintf(stderr, " - Fractional difference in log likelihoods = %lf%%\n", Lfrac);

  XLALFree(cdata);
  gsl_vector_complex_free(cmodelfull);
  gsl_vector_complex_free(cmodelreduced);
  gsl_matrix_complex_free(cmmw);
  gsl_vector_complex_free(cdmw);
  LALInferenceRemoveREALROQInterpolant( interp );
  LALInferenceRemoveCOMPLEXROQInterpolant( cinterp );

  /* check log likelihood difference is within tolerance */
  if ( Lfrac > LTOL ) { return 1; }

  return 0;
}
static int
XLALSimIMRSpinEOBInitialConditionsPrec(
				   REAL8Vector * initConds,	/**<< OUTPUT, Initial dynamical variables */
				   const REAL8 mass1,	/**<< mass 1 */
				   const REAL8 mass2,	/**<< mass 2 */
				   const REAL8 fMin,	/**<< Initial frequency (given) */
				   const REAL8 inc,	/**<< Inclination */
				   const REAL8 spin1[],	/**<< Initial spin vector 1 */
				   const REAL8 spin2[],	/**<< Initial spin vector 2 */
				   SpinEOBParams * params	/**<< Spin EOB parameters */
)
{

#ifndef LAL_NDEBUG
	if (!initConds) {
		XLAL_ERROR(XLAL_EINVAL);
	}
#endif

	int	debugPK = 0; int printPK = 0;
  FILE* UNUSED out = NULL;

	if (printPK) {
		XLAL_PRINT_INFO("Inside the XLALSimIMRSpinEOBInitialConditionsPrec function!\n");
		XLAL_PRINT_INFO(
    "Inputs: m1 = %.16e, m2 = %.16e, fMin = %.16e, inclination = %.16e\n",
      mass1, mass2, (double)fMin, (double)inc);
		XLAL_PRINT_INFO("Inputs: mSpin1 = {%.16e, %.16e, %.16e}\n",
      spin1[0], spin1[1], spin1[2]);
		XLAL_PRINT_INFO("Inputs: mSpin2 = {%.16e, %.16e, %.16e}\n",
      spin2[0], spin2[1], spin2[2]);
		fflush(NULL);
	}
	static const int UNUSED lMax = 8;

	int		i;

	/* Variable to keep track of whether the user requested the tortoise */
	int		tmpTortoise;

	UINT4		SpinAlignedEOBversion;

	REAL8		mTotal;
	REAL8		eta;
	REAL8		omega   , v0;	/* Initial velocity and angular
					 * frequency */

	REAL8		ham;	/* Hamiltonian */

	REAL8		LnHat    [3];	/* Initial orientation of angular
					 * momentum */
	REAL8		rHat     [3];	/* Initial orientation of radial
					 * vector */
	REAL8		vHat     [3];	/* Initial orientation of velocity
					 * vector */
	REAL8		Lhat     [3];	/* Direction of relativistic ang mom */
	REAL8		qHat     [3];
	REAL8		pHat     [3];

	/* q and p vectors in Cartesian and spherical coords */
	REAL8		qCart    [3], pCart[3];
	REAL8		qSph     [3], pSph[3];

	/* We will need to manipulate the spin vectors */
	/* We will use temporary vectors to do this */
	REAL8		tmpS1    [3];
	REAL8		tmpS2    [3];
	REAL8		tmpS1Norm[3];
	REAL8		tmpS2Norm[3];

	REAL8Vector	qCartVec, pCartVec;
	REAL8Vector	s1Vec, s2Vec, s1VecNorm, s2VecNorm;
	REAL8Vector	sKerr, sStar;
	REAL8		sKerrData[3], sStarData[3];
	REAL8		a = 0.;
	//, chiS, chiA;
	//REAL8 chi1, chi2;

	/*
	 * We will need a full values vector for calculating derivs of
	 * Hamiltonian
	 */
	REAL8		sphValues[12];
	REAL8		cartValues[12];

	/* Matrices for rotating to the new basis set. */
	/* It is more convenient to calculate the ICs in a simpler basis */
	gsl_matrix     *rotMatrix = NULL;
	gsl_matrix     *invMatrix = NULL;
	gsl_matrix     *rotMatrix2 = NULL;
	gsl_matrix     *invMatrix2 = NULL;

	/* Root finding stuff for finding the spherical orbit */
	SEOBRootParams	rootParams;
	const gsl_multiroot_fsolver_type *T = gsl_multiroot_fsolver_hybrid;
	gsl_multiroot_fsolver *rootSolver = NULL;

	gsl_multiroot_function rootFunction;
	gsl_vector     *initValues = NULL;
	gsl_vector     *finalValues = NULL;
	INT4 gslStatus;
        INT4 cntGslNoProgress = 0, MAXcntGslNoProgress = 5;
        //INT4 cntGslNoProgress = 0, MAXcntGslNoProgress = 50;
        REAL8 multFacGslNoProgress = 3./5.;
	//const int	maxIter = 2000;
	const int	maxIter = 10000;

	memset(&rootParams, 0, sizeof(rootParams));

	mTotal = mass1 + mass2;
	eta = mass1 * mass2 / (mTotal * mTotal);
	memcpy(tmpS1, spin1, sizeof(tmpS1));
	memcpy(tmpS2, spin2, sizeof(tmpS2));
	memcpy(tmpS1Norm, spin1, sizeof(tmpS1Norm));
	memcpy(tmpS2Norm, spin2, sizeof(tmpS2Norm));
	for (i = 0; i < 3; i++) {
		tmpS1Norm[i] /= mTotal * mTotal;
		tmpS2Norm[i] /= mTotal * mTotal;
	}
	SpinAlignedEOBversion = params->seobCoeffs->SpinAlignedEOBversion;
	/* We compute the ICs for the non-tortoise p, and convert at the end */
	tmpTortoise = params->tortoise;
	params->tortoise = 0;

	EOBNonQCCoeffs *nqcCoeffs = NULL;
	nqcCoeffs = params->nqcCoeffs;

	/*
	 * STEP 1) Rotate to LNhat0 along z-axis and N0 along x-axis frame,
	 * where LNhat0 and N0 are initial normal to orbital plane and
	 * initial orbital separation;
	 */

	/* Set the initial orbital ang mom direction. Taken from STPN code */
	LnHat[0] = sin(inc);
	LnHat[1] = 0.;
	LnHat[2] = cos(inc);

	/*
	 * Set the radial direction - need to take care to avoid singularity
	 * if L is along z axis
	 */
	if (LnHat[2] > 0.9999) {
		rHat[0] = 1.;
		rHat[1] = rHat[2] = 0.;
	} else {
		REAL8		theta0 = atan(-LnHat[2] / LnHat[0]);	/* theta0 is between 0
									 * and Pi */
		rHat[0] = sin(theta0);
		rHat[1] = 0;
		rHat[2] = cos(theta0);
	}

	/* Now we can complete the triad */
	vHat[0] = CalculateCrossProductPrec(0, LnHat, rHat);
	vHat[1] = CalculateCrossProductPrec(1, LnHat, rHat);
	vHat[2] = CalculateCrossProductPrec(2, LnHat, rHat);

	NormalizeVectorPrec(vHat);

	/* Vectors BEFORE rotation */
	if (printPK) {
		for (i = 0; i < 3; i++)
			XLAL_PRINT_INFO(" LnHat[%d] = %.16e, rHat[%d] = %.16e, vHat[%d] = %.16e\n",
        i, LnHat[i], i, rHat[i], i, vHat[i]);

		XLAL_PRINT_INFO("\n\n");
		for (i = 0; i < 3; i++)
			XLAL_PRINT_INFO(" s1[%d] = %.16e, s2[%d] = %.16e\n", i, tmpS1[i], i, tmpS2[i]);
    fflush(NULL);
	}

	/* Allocate and compute the rotation matrices */
	XLAL_CALLGSL(rotMatrix = gsl_matrix_alloc(3, 3));
	XLAL_CALLGSL(invMatrix = gsl_matrix_alloc(3, 3));
	if (!rotMatrix || !invMatrix) {
		if (rotMatrix)
			gsl_matrix_free(rotMatrix);
		if (invMatrix)
			gsl_matrix_free(invMatrix);
		XLAL_ERROR(XLAL_ENOMEM);
	}
	if (CalculateRotationMatrixPrec(rotMatrix, invMatrix, rHat, vHat, LnHat) == XLAL_FAILURE) {
		gsl_matrix_free(rotMatrix);
		gsl_matrix_free(invMatrix);
		XLAL_ERROR(XLAL_ENOMEM);
	}
	/* Rotate the orbital vectors and spins */
	ApplyRotationMatrixPrec(rotMatrix, rHat);
	ApplyRotationMatrixPrec(rotMatrix, vHat);
	ApplyRotationMatrixPrec(rotMatrix, LnHat);
	ApplyRotationMatrixPrec(rotMatrix, tmpS1);
	ApplyRotationMatrixPrec(rotMatrix, tmpS2);
	ApplyRotationMatrixPrec(rotMatrix, tmpS1Norm);
	ApplyRotationMatrixPrec(rotMatrix, tmpS2Norm);

	/* See if Vectors have been rotated fine */
	if (printPK) {
		XLAL_PRINT_INFO("\nAfter applying rotation matrix:\n\n");
		for (i = 0; i < 3; i++)
			XLAL_PRINT_INFO(" LnHat[%d] = %.16e, rHat[%d] = %.16e, vHat[%d] = %.16e\n",
                i, LnHat[i], i, rHat[i], i, vHat[i]);

		XLAL_PRINT_INFO("\n");
		for (i = 0; i < 3; i++)
			XLAL_PRINT_INFO(" s1[%d] = %.16e, s2[%d] = %.16e\n", i, tmpS1[i], i, tmpS2[i]);

    fflush(NULL);
	}
	/*
	 * STEP 2) After rotation in STEP 1, in spherical coordinates, phi0
	 * and theta0 are given directly in Eq. (4.7), r0, pr0, ptheta0 and
	 * pphi0 are obtained by solving Eqs. (4.8) and (4.9) (using
	 * gsl_multiroot_fsolver). At this step, we find initial conditions
	 * for a spherical orbit without radiation reaction.
	 */

  /* Initialise the gsl stuff */
	XLAL_CALLGSL(rootSolver = gsl_multiroot_fsolver_alloc(T, 3));
	if (!rootSolver) {
		gsl_matrix_free(rotMatrix);
		gsl_matrix_free(invMatrix);
		XLAL_ERROR(XLAL_ENOMEM);
	}
	XLAL_CALLGSL(initValues = gsl_vector_calloc(3));
	if (!initValues) {
		gsl_multiroot_fsolver_free(rootSolver);
		gsl_matrix_free(rotMatrix);
		gsl_matrix_free(invMatrix);
		XLAL_ERROR(XLAL_ENOMEM);
	}

	rootFunction.f = XLALFindSphericalOrbitPrec;
	rootFunction.n = 3;
	rootFunction.params = &rootParams;

	/* Calculate the initial velocity from the given initial frequency */
	omega = LAL_PI * mTotal * LAL_MTSUN_SI * fMin;
	v0 = cbrt(omega);

	/* Given this, we can start to calculate the initial conditions */
	/* for spherical coords in the new basis */
	rootParams.omega = omega;
	rootParams.params = params;

	/* To start with, we will just assign Newtonian-ish ICs to the system */
	rootParams.values[0] = scale1 * 1. / (v0 * v0);	/* Initial r */
	rootParams.values[4] = scale2 * v0;	            /* Initial p */
	rootParams.values[5] = scale3 * 1e-3;
	//PK
  memcpy(rootParams.values + 6, tmpS1, sizeof(tmpS1));
	memcpy(rootParams.values + 9, tmpS2, sizeof(tmpS2));

	if (printPK) {
    XLAL_PRINT_INFO("ICs guess: x = %.16e, py = %.16e, pz = %.16e\n",
      rootParams.values[0]/scale1, rootParams.values[4]/scale2,
      rootParams.values[5]/scale3);
    fflush(NULL);
  }

	gsl_vector_set(initValues, 0, rootParams.values[0]);
	gsl_vector_set(initValues, 1, rootParams.values[4]);
  gsl_vector_set(initValues, 2, rootParams.values[5]);

	gsl_multiroot_fsolver_set(rootSolver, &rootFunction, initValues);

	/* We are now ready to iterate to find the solution */
	i = 0;

  if(debugPK){ out = fopen("ICIterations.dat", "w"); }
	do {
		XLAL_CALLGSL(gslStatus = gsl_multiroot_fsolver_iterate(rootSolver));
		if (debugPK) {
      fprintf( out, "%d\t", i );

      /* Write to file */
      fprintf( out, "%.16e\t%.16e\t%.16e\t",
        rootParams.values[0]/scale1, rootParams.values[4]/scale2,
        rootParams.values[5]/scale3 );

      /* Residual Function values whose roots we are trying to find */
      finalValues = gsl_multiroot_fsolver_f(rootSolver);

      /* Write to file */
      fprintf( out, "%.16e\t%.16e\t%.16e\t",
        gsl_vector_get(finalValues, 0),
        gsl_vector_get(finalValues, 1),
        gsl_vector_get(finalValues, 2) );

      /* Step sizes in each of function variables */
      finalValues = gsl_multiroot_fsolver_dx(rootSolver);

      /* Write to file */
      fprintf( out, "%.16e\t%.16e\t%.16e\t%d\n",
        gsl_vector_get(finalValues, 0)/scale1,
        gsl_vector_get(finalValues, 1)/scale2,
        gsl_vector_get(finalValues, 2)/scale3,
        gslStatus );
		}

    if (gslStatus == GSL_ENOPROG || gslStatus == GSL_ENOPROGJ) {
      XLAL_PRINT_INFO(
        "\n NO PROGRESS being made by Spherical orbit root solver\n");

      /* Print Residual Function values whose roots we are trying to find */
      finalValues = gsl_multiroot_fsolver_f(rootSolver);
      XLAL_PRINT_INFO("Function value here given by the following:\n");
      XLAL_PRINT_INFO(" F1 = %.16e, F2 = %.16e, F3 = %.16e\n",
          gsl_vector_get(finalValues, 0),
		       gsl_vector_get(finalValues, 1), gsl_vector_get(finalValues, 2));

      /* Print Step sizes in each of function variables */
      finalValues = gsl_multiroot_fsolver_dx(rootSolver);
//      XLAL_PRINT_INFO("Stepsizes in each dimension:\n");
//      XLAL_PRINT_INFO(" x = %.16e, py = %.16e, pz = %.16e\n",
//          gsl_vector_get(finalValues, 0)/scale1,
//	       gsl_vector_get(finalValues, 1)/scale2,
//          gsl_vector_get(finalValues, 2)/scale3);

      /* Only allow this flag to be caught MAXcntGslNoProgress no. of times */
      cntGslNoProgress += 1;
      if (cntGslNoProgress >= MAXcntGslNoProgress) {
        cntGslNoProgress = 0;

        if(multFacGslNoProgress < 1.){ multFacGslNoProgress *= 1.02; }
        else{ multFacGslNoProgress /= 1.01; }

      } 
      /* Now that no progress is being made, we need to reset the initial guess
       * for the (r,pPhi, pTheta) and reset the integrator */
      rootParams.values[0] = scale1 * 1. / (v0 * v0);	/* Initial r */
      rootParams.values[4] = scale2 * v0;	            /* Initial p */
      if( cntGslNoProgress % 2 )
        rootParams.values[5] = scale3 * 1e-3 / multFacGslNoProgress;
      else
        rootParams.values[5] = scale3 * 1e-3 * multFacGslNoProgress;
      //PK
      memcpy(rootParams.values + 6, tmpS1, sizeof(tmpS1));
      memcpy(rootParams.values + 9, tmpS2, sizeof(tmpS2));

      if (printPK) {
        XLAL_PRINT_INFO("New ICs guess: x = %.16e, py = %.16e, pz = %.16e\n",
                rootParams.values[0]/scale1, rootParams.values[4]/scale2,
                rootParams.values[5]/scale3);
        fflush(NULL);
      }

      gsl_vector_set(initValues, 0, rootParams.values[0]);
      gsl_vector_set(initValues, 1, rootParams.values[4]);
      gsl_vector_set(initValues, 2, rootParams.values[5]);
      gsl_multiroot_fsolver_set(rootSolver, &rootFunction, initValues);
    }
    else if (gslStatus == GSL_EBADFUNC) {
      XLALPrintError(
      "Inf or Nan encountered in evaluluation of spherical orbit Eqn\n");
			gsl_multiroot_fsolver_free(rootSolver);
			gsl_vector_free(initValues);
			gsl_matrix_free(rotMatrix);
			gsl_matrix_free(invMatrix);
			XLAL_ERROR(XLAL_EDOM);
    }
		else if (gslStatus != GSL_SUCCESS) {
			XLALPrintError("Error in GSL iteration function!\n");
			gsl_multiroot_fsolver_free(rootSolver);
			gsl_vector_free(initValues);
			gsl_matrix_free(rotMatrix);
			gsl_matrix_free(invMatrix);
			XLAL_ERROR(XLAL_EDOM);
		}

    /* different ways to test convergence of the method */
		XLAL_CALLGSL(gslStatus = gsl_multiroot_test_residual(rootSolver->f, 1.0e-8));
    /*XLAL_CALLGSL(gslStatus= gsl_multiroot_test_delta(
          gsl_multiroot_fsolver_dx(rootSolver),
          gsl_multiroot_fsolver_root(rootSolver),
          1.e-8, 1.e-5));*/
		i++;
	}
	while (gslStatus == GSL_CONTINUE && i <= maxIter);

  if(debugPK) { fflush(NULL); fclose(out); }

	if (i > maxIter && gslStatus != GSL_SUCCESS) {
		gsl_multiroot_fsolver_free(rootSolver);
		gsl_vector_free(initValues);
		gsl_matrix_free(rotMatrix);
		gsl_matrix_free(invMatrix);
		//XLAL_ERROR(XLAL_EMAXITER);
		XLAL_ERROR(XLAL_EDOM);
	}
	finalValues = gsl_multiroot_fsolver_root(rootSolver);

	if (printPK) {
		XLAL_PRINT_INFO("Spherical orbit conditions here given by the following:\n");
		XLAL_PRINT_INFO(" x = %.16e, py = %.16e, pz = %.16e\n",
           gsl_vector_get(finalValues, 0)/scale1,
		       gsl_vector_get(finalValues, 1)/scale2,
           gsl_vector_get(finalValues, 2)/scale3);
	}
	memset(qCart, 0, sizeof(qCart));
	memset(pCart, 0, sizeof(pCart));

	qCart[0] = gsl_vector_get(finalValues, 0)/scale1;
	pCart[1] = gsl_vector_get(finalValues, 1)/scale2;
	pCart[2] = gsl_vector_get(finalValues, 2)/scale3;


	/* Free the GSL root finder, since we're done with it */
	gsl_multiroot_fsolver_free(rootSolver);
	gsl_vector_free(initValues);


	/*
	 * STEP 3) Rotate to L0 along z-axis and N0 along x-axis frame, where
	 * L0 is the initial orbital angular momentum and L0 is calculated
	 * using initial position and linear momentum obtained in STEP 2.
	 */

	/* Now we can calculate the relativistic L and rotate to a new basis */
	memcpy(qHat, qCart, sizeof(qCart));
	memcpy(pHat, pCart, sizeof(pCart));

	NormalizeVectorPrec(qHat);
	NormalizeVectorPrec(pHat);

	Lhat[0] = CalculateCrossProductPrec(0, qHat, pHat);
	Lhat[1] = CalculateCrossProductPrec(1, qHat, pHat);
	Lhat[2] = CalculateCrossProductPrec(2, qHat, pHat);

	NormalizeVectorPrec(Lhat);

	XLAL_CALLGSL(rotMatrix2 = gsl_matrix_alloc(3, 3));
	XLAL_CALLGSL(invMatrix2 = gsl_matrix_alloc(3, 3));

	if (CalculateRotationMatrixPrec(rotMatrix2, invMatrix2, qHat, pHat, Lhat) == XLAL_FAILURE) {
		gsl_matrix_free(rotMatrix);
		gsl_matrix_free(invMatrix);
		XLAL_ERROR(XLAL_ENOMEM);
	}
	ApplyRotationMatrixPrec(rotMatrix2, rHat);
	ApplyRotationMatrixPrec(rotMatrix2, vHat);
	ApplyRotationMatrixPrec(rotMatrix2, LnHat);
	ApplyRotationMatrixPrec(rotMatrix2, tmpS1);
	ApplyRotationMatrixPrec(rotMatrix2, tmpS2);
	ApplyRotationMatrixPrec(rotMatrix2, tmpS1Norm);
	ApplyRotationMatrixPrec(rotMatrix2, tmpS2Norm);
	ApplyRotationMatrixPrec(rotMatrix2, qCart);
	ApplyRotationMatrixPrec(rotMatrix2, pCart);

        gsl_matrix_free(rotMatrix);
        gsl_matrix_free(rotMatrix2);

        if (printPK) {
		XLAL_PRINT_INFO("qCart after rotation2 %3.10f %3.10f %3.10f\n", qCart[0], qCart[1], qCart[2]);
		XLAL_PRINT_INFO("pCart after rotation2 %3.10f %3.10f %3.10f\n", pCart[0], pCart[1], pCart[2]);
		XLAL_PRINT_INFO("S1 after rotation2 %3.10f %3.10f %3.10f\n", tmpS1Norm[0], tmpS1Norm[1], tmpS1Norm[2]);
		XLAL_PRINT_INFO("S2 after rotation2 %3.10f %3.10f %3.10f\n", tmpS2Norm[0], tmpS2Norm[1], tmpS2Norm[2]);
	}
	/*
	 * STEP 4) In the L0-N0 frame, we calculate (dE/dr)|sph using Eq.
	 * (4.14), then initial dr/dt using Eq. (4.10), and finally pr0 using
	 * Eq. (4.15).
	 */

	/* Now we can calculate the flux. Change to spherical co-ords */
	CartesianToSphericalPrec(qSph, pSph, qCart, pCart);
	memcpy(sphValues, qSph, sizeof(qSph));
	memcpy(sphValues + 3, pSph, sizeof(pSph));
	memcpy(sphValues + 6, tmpS1, sizeof(tmpS1));
	memcpy(sphValues + 9, tmpS2, sizeof(tmpS2));

	memcpy(cartValues, qCart, sizeof(qCart));
	memcpy(cartValues + 3, pCart, sizeof(pCart));
	memcpy(cartValues + 6, tmpS1, sizeof(tmpS1));
	memcpy(cartValues + 9, tmpS2, sizeof(tmpS2));

	REAL8		dHdpphi , d2Hdr2, d2Hdrdpphi;
	REAL8		rDot    , dHdpr, flux, dEdr;

	d2Hdr2 = XLALCalculateSphHamiltonianDeriv2Prec(0, 0, sphValues, params);
	d2Hdrdpphi = XLALCalculateSphHamiltonianDeriv2Prec(0, 5, sphValues, params);

	if (printPK)
		XLAL_PRINT_INFO("d2Hdr2 = %.16e, d2Hdrdpphi = %.16e\n", d2Hdr2, d2Hdrdpphi);

	/* New code to compute derivatives w.r.t. cartesian variables */

	REAL8		tmpDValues[14];
	int UNUSED	status;
	for (i = 0; i < 3; i++) {
		cartValues[i + 6] /= mTotal * mTotal;
		cartValues[i + 9] /= mTotal * mTotal;
	}
    UINT4 oldignoreflux = params->ignoreflux;
    params->ignoreflux = 1;
	status = XLALSpinPrecHcapNumericalDerivative(0, cartValues, tmpDValues, params);
    params->ignoreflux = oldignoreflux;
	for (i = 0; i < 3; i++) {
		cartValues[i + 6] *= mTotal * mTotal;
		cartValues[i + 9] *= mTotal * mTotal;
	}

	dHdpphi = tmpDValues[1] / sqrt(cartValues[0] * cartValues[0] + cartValues[1] * cartValues[1] + cartValues[2] * cartValues[2]);
	//XLALSpinPrecHcapNumDerivWRTParam(4, cartValues, params) / sphValues[0];

	dEdr = -dHdpphi * d2Hdr2 / d2Hdrdpphi;

	if (printPK)
		XLAL_PRINT_INFO("d2Hdr2 = %.16e d2Hdrdpphi = %.16e dHdpphi = %.16e\n",
            d2Hdr2, d2Hdrdpphi, dHdpphi);

	if (d2Hdr2 != 0.0) {
		/* We will need to calculate the Hamiltonian to get the flux */
		s1Vec.length = s2Vec.length = s1VecNorm.length = s2VecNorm.length = sKerr.length = sStar.length = 3;
		s1Vec.data = tmpS1;
		s2Vec.data = tmpS2;
		s1VecNorm.data = tmpS1Norm;
		s2VecNorm.data = tmpS2Norm;
		sKerr.data = sKerrData;
		sStar.data = sStarData;

		qCartVec.length = pCartVec.length = 3;
		qCartVec.data = qCart;
		pCartVec.data = pCart;

		//chi1 = tmpS1[0] * LnHat[0] + tmpS1[1] * LnHat[1] + tmpS1[2] * LnHat[2];
		//chi2 = tmpS2[0] * LnHat[0] + tmpS2[1] * LnHat[1] + tmpS2[2] * LnHat[2];

		//if (debugPK)
			//XLAL_PRINT_INFO("magS1 = %.16e, magS2 = %.16e\n", chi1, chi2);

		//chiS = 0.5 * (chi1 / (mass1 * mass1) + chi2 / (mass2 * mass2));
		//chiA = 0.5 * (chi1 / (mass1 * mass1) - chi2 / (mass2 * mass2));

		XLALSimIMRSpinEOBCalculateSigmaKerr(&sKerr, mass1, mass2, &s1Vec, &s2Vec);
		XLALSimIMRSpinEOBCalculateSigmaStar(&sStar, mass1, mass2, &s1Vec, &s2Vec);

		/*
		 * The a in the flux has been set to zero, but not in the
		 * Hamiltonian
		 */
		a = sqrt(sKerr.data[0] * sKerr.data[0] + sKerr.data[1] * sKerr.data[1] + sKerr.data[2] * sKerr.data[2]);
		//XLALSimIMREOBCalcSpinPrecFacWaveformCoefficients(params->eobParams->hCoeffs, mass1, mass2, eta, /* a */ 0.0, chiS, chiA);
		//XLALSimIMRCalculateSpinPrecEOBHCoeffs(params->seobCoeffs, eta, a);
		ham = XLALSimIMRSpinPrecEOBHamiltonian(eta, &qCartVec, &pCartVec, &s1VecNorm, &s2VecNorm, &sKerr, &sStar, params->tortoise, params->seobCoeffs);

		if (printPK)
			XLAL_PRINT_INFO("Stas: hamiltonian in ICs at this point is %.16e\n", ham);

		/* And now, finally, the flux */
		REAL8Vector	polarDynamics, cartDynamics;
		REAL8		polarData[4], cartData[12];

		polarDynamics.length = 4;
		polarDynamics.data = polarData;

		polarData[0] = qSph[0];
		polarData[1] = 0.;
		polarData[2] = pSph[0];
		polarData[3] = pSph[2];

		cartDynamics.length = 12;
		cartDynamics.data = cartData;

		memcpy(cartData, qCart, 3 * sizeof(REAL8));
		memcpy(cartData + 3, pCart, 3 * sizeof(REAL8));
		memcpy(cartData + 6, tmpS1Norm, 3 * sizeof(REAL8));
		memcpy(cartData + 9, tmpS2Norm, 3 * sizeof(REAL8));

		//XLAL_PRINT_INFO("Stas: starting FLux calculations\n");

		flux = XLALInspiralPrecSpinFactorizedFlux(&polarDynamics, &cartDynamics, nqcCoeffs, omega, params, ham, lMax, SpinAlignedEOBversion);
		/*
		 * flux  = XLALInspiralSpinFactorizedFlux( &polarDynamics,
		 * nqcCoeffs, omega, params, ham, lMax, SpinAlignedEOBversion
		 * );
		 */
		//XLAL_PRINT_INFO("Stas flux = %.16e \n", flux);
		//exit(0);
		flux = flux / eta;

		rDot = -flux / dEdr;
		if (debugPK) {
			XLAL_PRINT_INFO("Stas here I am 2  \n");
		}
		/*
		 * We now need dHdpr - we take it that it is safely linear up
		 * to a pr of 1.0e-3 PK: Ideally, the pr should be of the
		 * order of other momenta, in order for its contribution to
		 * the Hamiltonian to not get buried in the numerical noise
		 * in the numerically larger momenta components
		 */
		cartValues[3] = 1.0e-3;
		for (i = 0; i < 3; i++) {
			cartValues[i + 6] /= mTotal * mTotal;
			cartValues[i + 9] /= mTotal * mTotal;
		}
        oldignoreflux = params->ignoreflux;
        params->ignoreflux = 1;
        params->seobCoeffs->updateHCoeffs = 1;
		status = XLALSpinPrecHcapNumericalDerivative(0, cartValues, tmpDValues, params);
        params->ignoreflux = oldignoreflux;
		for (i = 0; i < 3; i++) {
			cartValues[i + 6] *= mTotal * mTotal;
			cartValues[i + 9] *= mTotal * mTotal;
		}
        REAL8		csi = sqrt(XLALSimIMRSpinPrecEOBHamiltonianDeltaT(params->seobCoeffs, qSph[0], eta, a)*XLALSimIMRSpinPrecEOBHamiltonianDeltaR(params->seobCoeffs, qSph[0], eta, a)) / (qSph[0] * qSph[0] + a * a);

		dHdpr = csi*tmpDValues[0];
		//XLALSpinPrecHcapNumDerivWRTParam(3, cartValues, params);

		if (debugPK) {
			XLAL_PRINT_INFO("Ingredients going into prDot:\n");
			XLAL_PRINT_INFO("flux = %.16e, dEdr = %.16e, dHdpr = %.16e, dHdpr/pr = %.16e\n", flux, dEdr, dHdpr, dHdpr / cartValues[3]);
		}
		/*
		 * We can now calculate what pr should be taking into account
		 * the flux
		 */
		pSph[0] = rDot / (dHdpr / cartValues[3]);
	} else {
		/*
		 * Since d2Hdr2 has evaluated to zero, we cannot do the
		 * above. Just set pr to zero
		 */
		//XLAL_PRINT_INFO("d2Hdr2 is zero!\n");
		pSph[0] = 0;
	}

	/* Now we are done - convert back to cartesian coordinates ) */
	SphericalToCartesianPrec(qCart, pCart, qSph, pSph);

	/*
	 * STEP 5) Rotate back to the original inertial frame by inverting
	 * the rotation of STEP 3 and then  inverting the rotation of STEP 1.
	 */

	/* Undo rotations to get back to the original basis */
	/* Second rotation */
	ApplyRotationMatrixPrec(invMatrix2, rHat);
	ApplyRotationMatrixPrec(invMatrix2, vHat);
	ApplyRotationMatrixPrec(invMatrix2, LnHat);
	ApplyRotationMatrixPrec(invMatrix2, tmpS1);
	ApplyRotationMatrixPrec(invMatrix2, tmpS2);
	ApplyRotationMatrixPrec(invMatrix2, tmpS1Norm);
	ApplyRotationMatrixPrec(invMatrix2, tmpS2Norm);
	ApplyRotationMatrixPrec(invMatrix2, qCart);
	ApplyRotationMatrixPrec(invMatrix2, pCart);

	/* First rotation */
	ApplyRotationMatrixPrec(invMatrix, rHat);
	ApplyRotationMatrixPrec(invMatrix, vHat);
	ApplyRotationMatrixPrec(invMatrix, LnHat);
	ApplyRotationMatrixPrec(invMatrix, tmpS1);
	ApplyRotationMatrixPrec(invMatrix, tmpS2);
	ApplyRotationMatrixPrec(invMatrix, tmpS1Norm);
	ApplyRotationMatrixPrec(invMatrix, tmpS2Norm);
	ApplyRotationMatrixPrec(invMatrix, qCart);
	ApplyRotationMatrixPrec(invMatrix, pCart);

        gsl_matrix_free(invMatrix);
        gsl_matrix_free(invMatrix2);

        /* If required, apply the tortoise transform */
	if (tmpTortoise) {
		REAL8		r = sqrt(qCart[0] * qCart[0] + qCart[1] * qCart[1] + qCart[2] * qCart[2]);
		REAL8		deltaR = XLALSimIMRSpinPrecEOBHamiltonianDeltaR(params->seobCoeffs, r, eta, a);
		REAL8		deltaT = XLALSimIMRSpinPrecEOBHamiltonianDeltaT(params->seobCoeffs, r, eta, a);
		REAL8		csi = sqrt(deltaT * deltaR) / (r * r + a * a);

		REAL8		pr = (qCart[0] * pCart[0] + qCart[1] * pCart[1] + qCart[2] * pCart[2]) / r;

		params->tortoise = tmpTortoise;

		if (debugPK) {
			XLAL_PRINT_INFO("Applying the tortoise to p (csi = %.26e)\n", csi);
			XLAL_PRINT_INFO("pCart = %3.10f %3.10f %3.10f\n", pCart[0], pCart[1], pCart[2]);
		}
		for (i = 0; i < 3; i++) {
			pCart[i] = pCart[i] + qCart[i] * pr * (csi - 1.) / r;
		}
	}


    /* Now copy the initial conditions back to the return vector */
	memcpy(initConds->data, qCart, sizeof(qCart));
	memcpy(initConds->data + 3, pCart, sizeof(pCart));
	memcpy(initConds->data + 6, tmpS1Norm, sizeof(tmpS1Norm));
	memcpy(initConds->data + 9, tmpS2Norm, sizeof(tmpS2Norm));

    for (i=0; i<12; i++) {
        if (fabs(initConds->data[i]) <=1.0e-15) {
            initConds->data[i] = 0.;
        }
    }

	if (debugPK) {
		XLAL_PRINT_INFO("THE FINAL INITIAL CONDITIONS:\n");
		XLAL_PRINT_INFO(" %.16e %.16e %.16e\n%.16e %.16e %.16e\n%.16e %.16e %.16e\n%.16e %.16e %.16e\n", initConds->data[0], initConds->data[1], initConds->data[2],
		       initConds->data[3], initConds->data[4], initConds->data[5], initConds->data[6], initConds->data[7], initConds->data[8],
		       initConds->data[9], initConds->data[10], initConds->data[11]);
	}
	return XLAL_SUCCESS;
}
/*
-------------------------------------------------------------------------
 * This function allocates and initialises the memory and parameters for
 * the workspace used for determining whether the two ellipsoids intersect.
 * If the shape matrices a and b are not null, these will be used in the
 * workspace; otherwise they will point to null, and will have to be
 * set by hand.
 * ------------------------------------------------------------------------*/
fContactWorkSpace * XLALInitFContactWorkSpace(
                       UINT4                          n,
                       const gsl_matrix              *a,
                       const gsl_matrix              *b,
                       const gsl_min_fminimizer_type *T,
                       REAL8                          conv
                                             )
{
  fContactWorkSpace *workSpace = NULL;

  /* Check the parameters passed in are sensible */
  if ( n <= 0 )
    XLAL_ERROR_NULL( XLAL_EINVAL );

  if ( conv <= 0.0 )
    XLAL_ERROR_NULL( XLAL_EINVAL );

  if ( !T )
    XLAL_ERROR_NULL( XLAL_EFAULT );

  /* The matrices a and b should be both supplied or both null */
  if ((!( a && b )) && ( a || b ))
        XLAL_ERROR_NULL( XLAL_EINVAL );

  /* Check the matrices conform to the expected sizes */
  if ( a )
  {
    if ( a->size1 != n || a->size2 != n || b->size1 != n || b->size2 != n )
    {
      XLAL_ERROR_NULL( XLAL_EBADLEN );
    }
  }


  /* Allocate the workspace */
  workSpace = (fContactWorkSpace *) LALCalloc( 1, sizeof(fContactWorkSpace) );
  if ( !workSpace )
    XLAL_ERROR_NULL( XLAL_ENOMEM );

  if ( a )
  {
    workSpace->invQ1         = a;
    workSpace->invQ2         = b;
  }

  XLAL_CALLGSL( workSpace->tmpA = gsl_matrix_calloc( n, n ) );
  XLAL_CALLGSL( workSpace->tmpB = gsl_matrix_calloc( n, n ) );
  XLAL_CALLGSL( workSpace->C    = gsl_matrix_calloc( n, n ) );
  XLAL_CALLGSL( workSpace->p1   = gsl_permutation_alloc( n ) );
  XLAL_CALLGSL( workSpace->tmpV = gsl_vector_calloc( n ) );
  XLAL_CALLGSL( workSpace->r_AB = gsl_vector_calloc( n ) );
  XLAL_CALLGSL( workSpace->s    = gsl_min_fminimizer_alloc( T ) );

  /* Check all of the above were allocated properly */
  if (!workSpace->tmpA || !workSpace->tmpB || !workSpace->C ||
      !workSpace->p1 || !workSpace->tmpV || !workSpace->r_AB
      || !workSpace->s )
  {
    XLALFreeFContactWorkSpace( workSpace );
    XLAL_ERROR_NULL( XLAL_ENOMEM );
  }

  /* Now set the rest of the parameters to the correct values */
  workSpace->T         = T;
  workSpace->convParam = conv;
  workSpace->n         = n;

  return workSpace;
}
/*
-------------------------------------------------------------------------
 * This function minimises fContact defined above using the
 * Brent method. It returns the minima with a negative sign (which then
 * becomes the maxima of the actual contact function. This can be compared
 * to 1 to check if two ellipsoids indeed overlap.
 * ------------------------------------------------------------------------*/
REAL8 XLALCheckOverlapOfEllipsoids (
        const gsl_vector   *ra,
        const gsl_vector   *rb,
        fContactWorkSpace  *workSpace )
{
    gsl_function        F;
    INT4                min_status;
    INT4                iter = 0, max_iter = 100;
    REAL8               m = 0.6180339887;
    REAL8               a = 0.0L, b = 1.0L;
    gsl_min_fminimizer  *s = workSpace->s;

    /* Sanity check on input arguments */
    if ( !ra || !rb || !workSpace )
      XLAL_ERROR_REAL8( XLAL_EFAULT );

    if ( ra->size != rb->size || ra->size != workSpace->n )
      XLAL_ERROR_REAL8( XLAL_EBADLEN);


    /* Set r_AB to be rb - ra */
    XLAL_CALLGSL( gsl_vector_memcpy( workSpace->r_AB, rb) );
    XLAL_CALLGSL( gsl_vector_sub (workSpace->r_AB, ra) );

    if ( gsl_vector_isnull( workSpace->r_AB ))
    {
      XLALPrintWarning("Position vectors ra and rb are identical.\n");
      return 0;
    }

    F.function = &fContact;
    F.params   = workSpace;

    XLAL_CALLGSL( min_status = gsl_min_fminimizer_set (s, &F, m, a, b) );
    if ( min_status != GSL_SUCCESS )
      XLAL_ERROR_REAL8( XLAL_EFUNC );

    do
    {
        iter++;
        XLAL_CALLGSL( min_status = gsl_min_fminimizer_iterate (s) );
        if (min_status != GSL_SUCCESS )
        {
          if (min_status == GSL_EBADFUNC)
            XLAL_ERROR_REAL8( XLAL_EFUNC | XLAL_EFPINVAL );
          else if (min_status == GSL_EZERODIV)
            XLAL_ERROR_REAL8( XLAL_EFUNC | XLAL_EFPDIV0 );
          else
            XLAL_ERROR_REAL8( XLAL_EFUNC );
        }

        m = gsl_min_fminimizer_x_minimum (s);
        a = gsl_min_fminimizer_x_lower (s);
        b = gsl_min_fminimizer_x_upper (s);

        XLAL_CALLGSL( min_status = gsl_min_test_interval (a, b, workSpace->convParam, 0.0) );
        if (min_status != GSL_CONTINUE && min_status != GSL_SUCCESS )
          XLAL_ERROR_REAL8( XLAL_EFUNC );
    }
    while (min_status == GSL_CONTINUE && iter < max_iter );
    /* End of minimization routine */

    /* Throw an error if max iterations would have been exceeded */
    if ( iter == max_iter && min_status == GSL_CONTINUE )
    {
      XLAL_ERROR_REAL8( XLAL_EMAXITER );
    }

    return ( -(s->f_minimum) );
}
REAL8 XLALMinimizeEThincaParameterOverTravelTime( REAL8 travelTime,
                                                  EThincaMinimizer *minimizer,
                                                  INT4   exttrig
                                                )
{
  REAL8 ethinca;


  /* If colocated detectors or known sky position, just return the e-thinca parameter */
  if (travelTime == 0.0 || exttrig )
  {
    ethinca = minimizeEThincaParameterOverTimeDiff( travelTime, minimizer );
    if ( XLAL_IS_REAL8_FAIL_NAN(ethinca) )
    {
      XLAL_ERROR_REAL8( XLAL_EFUNC );
    }
    return ethinca;
  }
  else
  {
    gsl_function        F;
    INT4                min_status;
    INT4                iter = 0;
    const INT4          max_iter = 100;
    REAL8               epsilon = 1.0 / 16384.0;
    REAL8               m = 0.0;
    REAL8               a = - travelTime, b = travelTime; /* Upper and lower bounds */
    REAL8               minEThinca, maxEThinca;
    REAL8               midEThinca;
    gsl_min_fminimizer  *s = gsl_min_fminimizer_alloc( minimizer->workSpace->T );

    if ( !s )
    {
      XLAL_ERROR_REAL8( XLAL_ENOMEM );
    }


    F.function = &minimizeEThincaParameterOverTimeDiff;
    F.params   = minimizer;

    /* Calculate e-thinca parameter at start, end and mid points */
    minEThinca = minimizeEThincaParameterOverTimeDiff( a, minimizer );
    maxEThinca = minimizeEThincaParameterOverTimeDiff( b, minimizer );
    midEThinca = minimizeEThincaParameterOverTimeDiff( m, minimizer );
    if ( XLAL_IS_REAL8_FAIL_NAN(minEThinca) || XLAL_IS_REAL8_FAIL_NAN(maxEThinca)
         || XLAL_IS_REAL8_FAIL_NAN(midEThinca) )
    {
      gsl_min_fminimizer_free( s );
      XLAL_ERROR_REAL8( XLAL_EFUNC );
    }

    /* Check we have contained a minimum. Otherwise take appropriate action */
    if ( midEThinca >= minEThinca || midEThinca >= maxEThinca )
    {
      REAL8 testEThinca; /* To contain the lowest end-point */
      if ( minEThinca < maxEThinca )
      {
        testEThinca = minEThinca;
        m           = a + 2.0 * epsilon;
      }
      else
      {
        testEThinca = maxEThinca;
        m = b - 2.0 * epsilon;
      }
      midEThinca = minimizeEThincaParameterOverTimeDiff( m, minimizer );
      if ( XLAL_IS_REAL8_FAIL_NAN(midEThinca) )
      {
        gsl_min_fminimizer_free( s );
        XLAL_ERROR_REAL8( XLAL_EFUNC );
      }

      /* If we still don't have the minimum return the lowest end-point */
      if ( midEThinca >= testEThinca )
      {
        gsl_min_fminimizer_free( s );
        return testEThinca;
      }
    }

    /* Set up the GSL minimizer */
    XLAL_CALLGSL( min_status = gsl_min_fminimizer_set_with_values(s, &F,
                       m, midEThinca, a, minEThinca, b, maxEThinca) );
    if ( min_status != GSL_SUCCESS )
    {
      gsl_min_fminimizer_free( s );
      XLAL_ERROR_REAL8( XLAL_EFUNC );
    }

    /* Loop to perform the minimization */
    do
    {
        iter++;
        XLAL_CALLGSL( min_status = gsl_min_fminimizer_iterate (s) );
        if (min_status != GSL_SUCCESS )
        {
            gsl_min_fminimizer_free( s );
            XLAL_ERROR_REAL8( XLAL_EFUNC );
        }

        m = gsl_min_fminimizer_x_minimum (s);
        a = gsl_min_fminimizer_x_lower (s);
        b = gsl_min_fminimizer_x_upper (s);

        XLAL_CALLGSL( min_status = gsl_min_test_interval (a, b, epsilon, 0.0) );
        if (min_status != GSL_CONTINUE && min_status != GSL_SUCCESS )
        {
          gsl_min_fminimizer_free( s );
          XLAL_ERROR_REAL8( XLAL_EFUNC );
        }
    }
    while ( min_status == GSL_CONTINUE && iter < max_iter );
    /* End of minimization routine */

    /* Throw an error if max iterations would have been exceeded */
    if ( iter == max_iter && min_status == GSL_CONTINUE )
    {
      gsl_min_fminimizer_free( s );
      XLAL_ERROR_REAL8( XLAL_EMAXITER );
    }

    /* Get the minimum e-thinca param, and free memory for minimizer */
    ethinca = gsl_min_fminimizer_f_minimum( s );
    gsl_min_fminimizer_free( s );
    XLALPrintInfo( "%s: Number of iterations = %d\n", __func__, iter);
  }

  /* Return the required e-thinca value */
  return ethinca;
}
/**
 * Function to calculate R.H.S. of the ODEs, given dyanmical variables,
 * their derivatives and EOB parameters. Since SEOBNRv1 spin Hamiltonian
 * was implemented for Cartesean coordinates while dynamical evolution was
 * implemented in polar coordinates, we need to perform a transform.
 * This is done in a particular transform in which
 * x = r, y = z = 0, px = pr, py = pphi/r, pz = 0, and
 * omega = v/r = (dy/dt)/r = (dH/dpy)/r, dr/dt = dx/dt = dH/dpx, etc.
 */
static int XLALSpinAlignedHcapDerivative(
                  double UNUSED t,          /**< UNUSED */
                  const REAL8   values[],   /**< dynamical varables */
                  REAL8         dvalues[],  /**< time derivative of dynamical variables */
                  void         *funcParams  /**< EOB parameters */
                  )
{

  static const REAL8 STEP_SIZE = 1.0e-4;

  static const INT4 lMax = 8;

  HcapDerivParams params;

  /* Since we take numerical derivatives wrt dynamical variables */
  /* but we want them wrt time, we use this temporary vector in  */
  /* the conversion */
  REAL8           tmpDValues[6];

  /* Cartesian values for calculating the Hamiltonian */
  REAL8           cartValues[6];

  REAL8           H; //Hamiltonian
  REAL8           flux;

  gsl_function F;
  INT4         gslStatus;
  UINT4 SpinAlignedEOBversion;
  UINT4 i;

  REAL8Vector rVec, pVec;
  REAL8 rData[3], pData[3];

  /* We need r, phi, pr, pPhi to calculate the flux */
  REAL8       r;
  REAL8Vector polarDynamics;
  REAL8       polData[4];

  REAL8 mass1, mass2, eta;

  /* Spins */
  REAL8Vector *s1Vec = NULL;
  REAL8Vector *s2Vec = NULL;
  REAL8Vector *sKerr = NULL;
  REAL8Vector *sStar = NULL;

  REAL8 a;

  REAL8 omega;

  /* EOB potential functions */
  REAL8 DeltaT, DeltaR;
  REAL8 csi;

  /* The error in a derivative as measured by GSL */
  REAL8 absErr;

  /* Declare NQC coefficients */
  EOBNonQCCoeffs *nqcCoeffs = NULL;

  /* Set up pointers for GSL */ 
  params.values  = cartValues;
  params.params  = (SpinEOBParams *)funcParams;
  nqcCoeffs = params.params->nqcCoeffs;

  s1Vec = params.params->s1Vec;
  s2Vec = params.params->s2Vec;
  sKerr = params.params->sigmaKerr;
  sStar = params.params->sigmaStar;

  F.function = &GSLSpinAlignedHamiltonianWrapper;
  F.params   = &params;

  mass1 = params.params->eobParams->m1;
  mass2 = params.params->eobParams->m2;
  eta   = params.params->eobParams->eta;

  SpinAlignedEOBversion = params.params->seobCoeffs->SpinAlignedEOBversion;

  r = values[0];

  /* Since this is spin aligned, I make the assumption */
  /* that the spin vector is along the z-axis.         */
  a  = sKerr->data[2];

  /* Calculate the potential functions and the tortoise coordinate factor csi,
     given by Eq. 28 of Pan et al. PRD 81, 084041 (2010) */
  DeltaT = XLALSimIMRSpinEOBHamiltonianDeltaT( params.params->seobCoeffs, r, eta, a );
  DeltaR = XLALSimIMRSpinEOBHamiltonianDeltaR( params.params->seobCoeffs, r, eta, a );
  csi    = sqrt( DeltaT * DeltaR ) / (r*r + a*a);
  //printf("DeltaT = %.16e, DeltaR = %.16e, a = %.16e\n",DeltaT,DeltaR,a);
  //printf( "csi in derivatives function = %.16e\n", csi );

  /* Populate the Cartesian values vector, using polar coordinate values */
  /* We can assume phi is zero wlog */
  memset( cartValues, 0, sizeof( cartValues ) );
  cartValues[0] = values[0];
  cartValues[3] = values[2];
  cartValues[4] = values[3] / values[0];

  /* Now calculate derivatives w.r.t. each Cartesian variable */
  for ( i = 0; i < 6; i++ )
  {
    params.varyParam = i;
    XLAL_CALLGSL( gslStatus = gsl_deriv_central( &F, cartValues[i], 
                    STEP_SIZE, &tmpDValues[i], &absErr ) );

    if ( gslStatus != GSL_SUCCESS )
    {
      XLALPrintError( "XLAL Error - %s: Failure in GSL function\n", __func__ );
      XLAL_ERROR( XLAL_EFUNC );
    }
  }

  /* Calculate the Cartesian vectors rVec and pVec */
  polarDynamics.length = 4;
  polarDynamics.data   = polData;

  memcpy( polData, values, sizeof( polData ) );

  rVec.length = pVec.length = 3;
  rVec.data   = rData;
  pVec.data   = pData;

  memset( rData, 0, sizeof(rData) );
  memset( pData, 0, sizeof(pData) );

  rData[0] = values[0];
  pData[0] = values[2];
  pData[1] = values[3] / values[0];
  /* Calculate Hamiltonian using Cartesian vectors rVec and pVec */
  H =  XLALSimIMRSpinEOBHamiltonian( eta, &rVec, &pVec, s1Vec, s2Vec, sKerr, sStar, params.params->tortoise, params.params->seobCoeffs );

  //printf( "csi = %.16e, ham = %.16e ( tortoise = %d)\n", csi, H, params.params->tortoise );
  //exit(1);
  //if ( values[0] > 1.3 && values[0] < 3.9 ) printf( "r = %e\n", values[0] );
  //if ( values[0] > 1.3 && values[0] < 3.9 ) printf( "Hamiltonian = %e\n", H );
  H = H * (mass1 + mass2);


  /*if ( values[0] > 1.3 && values[0] < 3.9 ) printf( "Cartesian derivatives:\n%f %f %f %f %f %f\n",
      tmpDValues[3], tmpDValues[4], tmpDValues[5], -tmpDValues[0], -tmpDValues[1], -tmpDValues[2] );*/

  /* Now calculate omega, and hence the flux */
  omega = tmpDValues[4] / r;
  flux  = XLALInspiralSpinFactorizedFlux( &polarDynamics, nqcCoeffs, omega, params.params, H/(mass1+mass2), lMax, SpinAlignedEOBversion );

  /* Looking at the non-spinning model, I think we need to divide the flux by eta */
  flux = flux / eta;

  //printf( "Flux in derivatives function = %.16e\n", flux );

  /* Now we can calculate the final (spherical) derivatives */
  /* csi is needed because we use the tortoise co-ordinate */
  /* Right hand side of Eqs. 10a - 10d of Pan et al. PRD 84, 124052 (2011) */
  dvalues[0] = csi * tmpDValues[3];
  dvalues[1] = omega;
  /* Note: in this special coordinate setting, namely y = z = 0, dpr/dt = dpx/dt + dy/dt * py/r, where py = pphi/r */ 
  dvalues[2] = - tmpDValues[0] + tmpDValues[4] * values[3] / (r*r);
  dvalues[2] = dvalues[2] * csi - ( values[2] / values[3] ) * flux / omega;
  dvalues[3] = - flux / omega;

  //if ( values[0] > 1.3 && values[0] < 3.9 ) printf("Values:\n%f %f %f %f\n", values[0], values[1], values[2], values[3] );

  //if ( values[0] > 1.3 && values[0] < 3.9 ) printf("Derivatives:\n%f %f %f %f\n", dvalues[0], r*dvalues[1], dvalues[2], dvalues[3] );

  if ( isnan( dvalues[0] ) || isnan( dvalues[1] ) || isnan( dvalues[2] ) || isnan( dvalues[3] ) )
  {
    //printf( "Deriv is nan: %e %e %e %e\n", dvalues[0], dvalues[1], dvalues[2], dvalues[3] );
    return 1;
  }

  return XLAL_SUCCESS;
}
/**
 * Generates the ringdown wave associated with the given real
 * and imaginary parts of the inspiral waveform. The parameters of
 * the ringdown, such as amplitude and phase offsets, are determined
 * by solving the linear equations defined in the DCC document T1100433.
 * In the linear equations Ax=y,
 * A is a 16-by-16 matrix depending on QNM (complex) frequencies,
 * x is a 16-d vector of the 8 unknown complex QNM amplitudes,
 * y is a 16-d vector depending on inspiral-plunge waveforms and their derivatives near merger.
 */
static INT4 XLALSimIMREOBHybridRingdownWave(
  REAL8Vector          *rdwave1,   /**<< OUTPUT, Real part of ringdown waveform */
  REAL8Vector          *rdwave2,   /**<< OUTPUT, Imag part of ringdown waveform */
  const REAL8           dt,        /**<< Sampling interval */
  const REAL8           mass1,     /**<< First component mass (in Solar masses) */
  const REAL8           mass2,     /**<< Second component mass (in Solar masses) */
  REAL8VectorSequence  *inspwave1, /**<< Values and derivs of real part inspiral waveform */
  REAL8VectorSequence  *inspwave2, /**<< Values and derivs of imag part inspiral waveform */
  COMPLEX16Vector      *modefreqs, /**<< Complex freqs of ringdown (scaled by total mass) */
  REAL8Vector          *matchrange /**<< Times which determine the comb of ringdown attachment */
  )
{

  /* XLAL error handling */
  INT4 errcode = XLAL_SUCCESS;

  /* For checking GSL return codes */
  INT4 gslStatus;

  UINT4 i, j, k, nmodes = 8;

  /* Sampling rate from input */
  REAL8 t1, t2, t3, t4, t5, rt;
  gsl_matrix *coef;
  gsl_vector *hderivs;
  gsl_vector *x;
  gsl_permutation *p;
  REAL8Vector *modeamps;
  int s;
  REAL8 tj;
  REAL8 m;

  /* mass in geometric units */
  m  = (mass1 + mass2) * LAL_MTSUN_SI;
  t5 = (matchrange->data[0] - matchrange->data[1]) * m;
  rt = -t5 / 5.;

  t4 = t5 + rt;
  t3 = t4 + rt;
  t2 = t3 + rt;
  t1 = t2 + rt;
  
  if ( inspwave1->length != 3 || inspwave2->length != 3 ||
		modefreqs->length != nmodes )
  {
    XLAL_ERROR( XLAL_EBADLEN );
  }

  /* Solving the linear system for QNMs amplitude coefficients using gsl routine */
  /* Initiate matrices and supporting variables */
  XLAL_CALLGSL( coef = (gsl_matrix *) gsl_matrix_alloc(2 * nmodes, 2 * nmodes) );
  XLAL_CALLGSL( hderivs = (gsl_vector *) gsl_vector_alloc(2 * nmodes) );
  XLAL_CALLGSL( x = (gsl_vector *) gsl_vector_alloc(2 * nmodes) );
  XLAL_CALLGSL( p = (gsl_permutation *) gsl_permutation_alloc(2 * nmodes) );

  /* Check all matrices and variables were allocated */
  if ( !coef || !hderivs || !x || !p )
  {
    if (coef)    gsl_matrix_free(coef);
    if (hderivs) gsl_vector_free(hderivs);
    if (x)       gsl_vector_free(x);
    if (p)       gsl_permutation_free(p);

    XLAL_ERROR( XLAL_ENOMEM );
  }

  /* Define the linear system Ax=y */
  /* Matrix A (2*n by 2*n) has block symmetry. Define half of A here as "coef" */
  /* The half of A defined here corresponds to matrices M1 and -M2 in the DCC document T1100433 */ 
  /* Define y here as "hderivs" */
  for (i = 0; i < nmodes; ++i)
  {
	gsl_matrix_set(coef, 0, i, 1);
	gsl_matrix_set(coef, 1, i, - cimag(modefreqs->data[i]));
	gsl_matrix_set(coef, 2, i, exp(-cimag(modefreqs->data[i])*t1) * cos(creal(modefreqs->data[i])*t1));
	gsl_matrix_set(coef, 3, i, exp(-cimag(modefreqs->data[i])*t2) * cos(creal(modefreqs->data[i])*t2));
	gsl_matrix_set(coef, 4, i, exp(-cimag(modefreqs->data[i])*t3) * cos(creal(modefreqs->data[i])*t3));
	gsl_matrix_set(coef, 5, i, exp(-cimag(modefreqs->data[i])*t4) * cos(creal(modefreqs->data[i])*t4));
	gsl_matrix_set(coef, 6, i, exp(-cimag(modefreqs->data[i])*t5) * cos(creal(modefreqs->data[i])*t5));
	gsl_matrix_set(coef, 7, i, exp(-cimag(modefreqs->data[i])*t5) * 
				      (-cimag(modefreqs->data[i]) * cos(creal(modefreqs->data[i])*t5)
				       -creal(modefreqs->data[i]) * sin(creal(modefreqs->data[i])*t5)));
	gsl_matrix_set(coef, 8, i, 0);
	gsl_matrix_set(coef, 9, i, - creal(modefreqs->data[i]));
	gsl_matrix_set(coef, 10, i, -exp(-cimag(modefreqs->data[i])*t1) * sin(creal(modefreqs->data[i])*t1));
	gsl_matrix_set(coef, 11, i, -exp(-cimag(modefreqs->data[i])*t2) * sin(creal(modefreqs->data[i])*t2));
	gsl_matrix_set(coef, 12, i, -exp(-cimag(modefreqs->data[i])*t3) * sin(creal(modefreqs->data[i])*t3));
	gsl_matrix_set(coef, 13, i, -exp(-cimag(modefreqs->data[i])*t4) * sin(creal(modefreqs->data[i])*t4));
	gsl_matrix_set(coef, 14, i, -exp(-cimag(modefreqs->data[i])*t5) * sin(creal(modefreqs->data[i])*t5));
	gsl_matrix_set(coef, 15, i, exp(-cimag(modefreqs->data[i])*t5) * 
				      ( cimag(modefreqs->data[i]) * sin(creal(modefreqs->data[i])*t5)
				       -creal(modefreqs->data[i]) * cos(creal(modefreqs->data[i])*t5)));
  }
  for (i = 0; i < 2; ++i)
  {
	gsl_vector_set(hderivs, i, inspwave1->data[(i + 1) * inspwave1->vectorLength - 1]);
	gsl_vector_set(hderivs, i + nmodes, inspwave2->data[(i + 1) * inspwave2->vectorLength - 1]);
	gsl_vector_set(hderivs, i + 6, inspwave1->data[i * inspwave1->vectorLength]);
	gsl_vector_set(hderivs, i + 6 + nmodes, inspwave2->data[i * inspwave2->vectorLength]);
  }
  gsl_vector_set(hderivs, 2, inspwave1->data[4]);
  gsl_vector_set(hderivs, 2 + nmodes, inspwave2->data[4]);
  gsl_vector_set(hderivs, 3, inspwave1->data[3]);
  gsl_vector_set(hderivs, 3 + nmodes, inspwave2->data[3]);
  gsl_vector_set(hderivs, 4, inspwave1->data[2]);
  gsl_vector_set(hderivs, 4 + nmodes, inspwave2->data[2]);
  gsl_vector_set(hderivs, 5, inspwave1->data[1]);
  gsl_vector_set(hderivs, 5 + nmodes, inspwave2->data[1]);
  
  /* Complete the definition for the rest half of A */
  for (i = 0; i < nmodes; ++i)
  {
	for (k = 0; k < nmodes; ++k)
	{
	  gsl_matrix_set(coef, i, k + nmodes, - gsl_matrix_get(coef, i + nmodes, k));
	  gsl_matrix_set(coef, i + nmodes, k + nmodes, gsl_matrix_get(coef, i, k));
	}
  }

  #if 0
  /* print ringdown-matching linear system: coefficient matrix and RHS vector */
  printf("\nRingdown matching matrix:\n");
  for (i = 0; i < 16; ++i)
  {
    for (j = 0; j < 16; ++j)
    {
      printf("%.12e ",gsl_matrix_get(coef,i,j));
    }
    printf("\n");
  }
  printf("RHS:  ");
  for (i = 0; i < 16; ++i)
  {
    printf("%.12e   ",gsl_vector_get(hderivs,i));
  }
  printf("\n");
  #endif

  /* Call gsl LU decomposition to solve the linear system */
  XLAL_CALLGSL( gslStatus = gsl_linalg_LU_decomp(coef, p, &s) );
  if ( gslStatus == GSL_SUCCESS )
  {
    XLAL_CALLGSL( gslStatus = gsl_linalg_LU_solve(coef, p, hderivs, x) );
  }
  if ( gslStatus != GSL_SUCCESS )
  {
    gsl_matrix_free(coef);
    gsl_vector_free(hderivs);
    gsl_vector_free(x);
    gsl_permutation_free(p);
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* Putting solution to an XLAL vector */
  modeamps = XLALCreateREAL8Vector(2 * nmodes);

  if ( !modeamps )
  {
    gsl_matrix_free(coef);
    gsl_vector_free(hderivs);
    gsl_vector_free(x);
    gsl_permutation_free(p);
    XLAL_ERROR( XLAL_ENOMEM );
  }

  for (i = 0; i < nmodes; ++i)
  {
	modeamps->data[i] = gsl_vector_get(x, i);
	modeamps->data[i + nmodes] = gsl_vector_get(x, i + nmodes);
  }

  /* Free all gsl linear algebra objects */
  gsl_matrix_free(coef);
  gsl_vector_free(hderivs);
  gsl_vector_free(x);
  gsl_permutation_free(p);

  /* Build ring-down waveforms */

  REAL8 timeOffset = fmod( matchrange->data[1], dt/m) * dt;

  for (j = 0; j < rdwave1->length; ++j)
  {
	tj = j * dt - timeOffset;
	rdwave1->data[j] = 0;
	rdwave2->data[j] = 0;
	for (i = 0; i < nmodes; ++i)
	{
	  rdwave1->data[j] += exp(- tj * cimag(modefreqs->data[i]))
			* ( modeamps->data[i] * cos(tj * creal(modefreqs->data[i]))
			+   modeamps->data[i + nmodes] * sin(tj * creal(modefreqs->data[i])) );
	  rdwave2->data[j] += exp(- tj * cimag(modefreqs->data[i]))
			* (- modeamps->data[i] * sin(tj * creal(modefreqs->data[i]))
			+   modeamps->data[i + nmodes] * cos(tj * creal(modefreqs->data[i])) );
	}
  }

  XLALDestroyREAL8Vector(modeamps);
  return errcode;
}
/**
 * Function which calculates the value of the waveform, plus its
 * first and second derivatives, for the points which will be required
 * in the hybrid comb attachment of the ringdown.
 */
static INT4 XLALGenerateHybridWaveDerivatives (
	REAL8Vector	*rwave,      /**<< OUTPUT, values of the waveform at comb points */
	REAL8Vector	*dwave,      /**<< OUTPUT, 1st deriv of the waveform at comb points */
	REAL8Vector	*ddwave,     /**<< OUTPUT, 2nd deriv of the waveform at comb points */
        REAL8Vector	*timeVec,    /**<< Vector containing the time */
	REAL8Vector	*wave,       /**<< Last part of inspiral waveform */
	REAL8Vector	*matchrange, /**<< Times which determine the size of the comb */
        REAL8           dt,          /**<< Sample time step */
        REAL8           mass1,       /**<< First component mass (in Solar masses) */
        REAL8           mass2        /**<< Second component mass (in Solar masses) */
	)
{

  /* XLAL error handling */
  INT4 errcode = XLAL_SUCCESS;

  /* For checking GSL return codes */
  INT4 gslStatus;

  UINT4 j;
  UINT4 vecLength;
  REAL8 m;
  double *y;
  double ry, dy, dy2;
  double rt;
  double *tlist;
  gsl_interp_accel *acc;
  gsl_spline *spline;

  /* Total mass in geometric units */
  m  = (mass1 + mass2) * LAL_MTSUN_SI;

  tlist = (double *) LALMalloc(6 * sizeof(double));
  rt = (matchrange->data[1] - matchrange->data[0]) / 5.;
  tlist[0] = matchrange->data[0];
  tlist[1] = tlist[0] + rt;
  tlist[2] = tlist[1] + rt;
  tlist[3] = tlist[2] + rt;
  tlist[4] = tlist[3] + rt;
  tlist[5] = matchrange->data[1];

  /* Set the length of the interpolation vectors */
  vecLength = ( m * matchrange->data[2] / dt ) + 1;

  /* Getting interpolation and derivatives of the waveform using gsl spline routine */
  /* Initiate arrays and supporting variables for gsl */
  y = (double *) LALMalloc(vecLength * sizeof(double));

  if ( !y )
  {
    XLAL_ERROR( XLAL_ENOMEM );
  }

  for (j = 0; j < vecLength; ++j)
  {
	y[j] = wave->data[j];
  }


  XLAL_CALLGSL( acc = (gsl_interp_accel*) gsl_interp_accel_alloc() );
  XLAL_CALLGSL( spline = (gsl_spline*) gsl_spline_alloc(gsl_interp_cspline, vecLength) );
  if ( !acc || !spline )
  {
    if ( acc )    gsl_interp_accel_free(acc);
    if ( spline ) gsl_spline_free(spline);
    LALFree( y );
    XLAL_ERROR( XLAL_ENOMEM );
  }

  /* Gall gsl spline interpolation */
  gslStatus = gsl_spline_init(spline, timeVec->data, y, vecLength);
  if ( gslStatus != GSL_SUCCESS )
  { 
    gsl_spline_free(spline);
    gsl_interp_accel_free(acc);
    LALFree( y );
    XLAL_ERROR( XLAL_EFUNC );
  }

  /* Getting first and second order time derivatives from gsl interpolations */
  for (j = 0; j < 6; ++j)
  {
    gslStatus = gsl_spline_eval_e( spline, tlist[j], acc, &ry );
    if ( gslStatus == GSL_SUCCESS )
    {
      gslStatus = gsl_spline_eval_deriv_e(spline, tlist[j], acc, &dy );
      gslStatus = gsl_spline_eval_deriv2_e(spline, tlist[j], acc, &dy2 );
    }
    if (gslStatus != GSL_SUCCESS )
    {
      gsl_spline_free(spline);
      gsl_interp_accel_free(acc);
      LALFree( y );
      XLAL_ERROR( XLAL_EFUNC );
    }
    rwave->data[j]  = (REAL8)(ry);
    dwave->data[j]  = (REAL8)(dy/m);
    ddwave->data[j] = (REAL8)(dy2/m/m);

  }
  
  /* Free gsl variables */
  gsl_spline_free(spline);
  gsl_interp_accel_free(acc);
  LALFree( tlist );
  LALFree(y);

  return errcode;
}