/**
 * \author Creighton, T. D.
 *
 * \brief Computes a continuous waveform with frequency drift and Doppler
 * modulation from a parabolic orbital trajectory.
 *
 * This function computes a quaiperiodic waveform using the spindown and
 * orbital parameters in <tt>*params</tt>, storing the result in
 * <tt>*output</tt>.
 *
 * In the <tt>*params</tt> structure, the routine uses all the "input"
 * fields specified in \ref GenerateSpinOrbitCW_h, and sets all of the
 * "output" fields.  If <tt>params-\>f</tt>=\c NULL, no spindown
 * modulation is performed.  If <tt>params-\>oneMinusEcc</tt>\f$\neq0\f$, or if
 * <tt>params-\>rPeriNorm</tt>\f$\times\f$<tt>params-\>angularSpeed</tt>\f$\geq1\f$
 * (faster-than-light speed at periapsis), an error is returned.
 *
 * In the <tt>*output</tt> structure, the field <tt>output-\>h</tt> is
 * ignored, but all other pointer fields must be set to \c NULL.  The
 * function will create and allocate space for <tt>output-\>a</tt>,
 * <tt>output-\>f</tt>, and <tt>output-\>phi</tt> as necessary.  The
 * <tt>output-\>shift</tt> field will remain set to \c NULL.
 *
 * ### Algorithm ###
 *
 * For parabolic orbits, we combine \eqref{eq_spinorbit-tr},
 * \eqref{eq_spinorbit-t}, and \eqref{eq_spinorbit-upsilon} to get \f$t_r\f$
 * directly as a function of \f$E\f$:
 * \f{equation}{
 * \label{eq_cubic-e}
 * t_r = t_p + \frac{r_p\sin i}{c} \left[ \cos\omega +
 * \left(\frac{1}{v_p} + \cos\omega\right)E -
 * \frac{\sin\omega}{4}E^2 + \frac{1}{12v_p}E^3\right] \;,
 * \f}
 * where \f$v_p=r_p\dot{\upsilon}_p\sin i/c\f$ is a normalized velocity at
 * periapsis.  Following the prescription for the general analytic
 * solution to the real cubic equation, we substitute
 * \f$E=x+3v_p\sin\omega\f$ to obtain:
 * \f{equation}{
 * \label{eq_cubic-x}
 * x^3 + px = q \;,
 * \f}
 * where:
 * \f{eqnarray}{
 * \label{eq_cubic-p}
 * p & = & 12 + 12v_p\cos\omega - 3v_p^2\sin^2\omega \;, \\
 * \label{eq_cubic-q}
 * q & = & 12v_p^2\sin\omega\cos\omega - 24v_p\sin\omega +
 * 2v_p^3\sin^3\omega + 12\dot{\upsilon}_p(t_r-t_p) \;.
 * \f}
 * We note that \f$p>0\f$ is guaranteed as long as \f$v_p<1\f$, so the right-hand
 * side of \eqref{eq_cubic-x} is monotonic in \f$x\f$ and has exactly one
 * root.  However, \f$p\rightarrow0\f$ in the limit \f$v_p\rightarrow1\f$ and
 * \f$\omega=\pi\f$.  This may cause some loss of precision in subsequent
 * calculations.  But \f$v_p\sim1\f$ means that our solution will be
 * inaccurate anyway because we ignore significant relativistic effects.
 *
 * Since \f$p>0\f$, we can substitute \f$x=y\sqrt{3/4p}\f$ to obtain:
 * \f{equation}{
 * 4y^3 + 3y = \frac{q}{2}\left(\frac{3}{p}\right)^{3/2} \equiv C \;.
 * \f}
 * Using the triple-angle hyperbolic identity
 * \f$\sinh(3\theta)=4\sinh^3\theta+3\sinh\theta\f$, we have
 * \f$y=\sinh\left(\frac{1}{3}\sinh^{-1}C\right)\f$.  The solution to the
 * original cubic equation is then:
 * \f{equation}{
 * E = 3v_p\sin\omega + 2\sqrt{\frac{p}{3}}
 * \sinh\left(\frac{1}{3}\sinh^{-1}C\right) \;.
 * \f}
 * To ease the calculation of \f$E\f$, we precompute the constant part
 * \f$E_0=3v_p\sin\omega\f$ and the coefficient \f$\Delta E=2\sqrt{p/3}\f$.
 * Similarly for \f$C\f$, we precompute a constant piece \f$C_0\f$ evaluated at
 * the epoch of the output time series, and a stepsize coefficient
 * \f$\Delta C=6(p/3)^{3/2}\dot{\upsilon}_p\Delta t\f$, where \f$\Delta t\f$ is
 * the step size in the (output) time series in \f$t_r\f$.  Thus at any
 * timestep \f$i\f$, we obtain \f$C\f$ and hence \f$E\f$ via:
 * \f{eqnarray}{
 * C & = & C_0 + i\Delta C \;,\\
 * E & = & E_0 + \Delta E\times\left\{\begin{array}{l@{\qquad}c}
 * \sinh\left[\frac{1}{3}\ln\left(
 * C + \sqrt{C^2+1} \right) \right]\;, & C\geq0 \;,\\ \\
 * \sinh\left[-\frac{1}{3}\ln\left(
 * -C + \sqrt{C^2+1} \right) \right]\;, & C\leq0 \;,\\
 * \end{array}\right.
 * \f}
 * where we have explicitly written \f$\sinh^{-1}\f$ in terms of functions in
 * \c math.h.  Once \f$E\f$ is found, we can compute
 * \f$t=E(12+E^2)/(12\dot{\upsilon}_p)\f$ (where again \f$1/12\dot{\upsilon}_p\f$
 * can be precomputed), and hence \f$f\f$ and \f$\phi\f$ via
 * \eqref{eq_taylorcw-freq} and \eqref{eq_taylorcw-phi}.  The
 * frequency \f$f\f$ must then be divided by the Doppler factor:
 * \f[
 * 1 + \frac{\dot{R}}{c} = 1 + \frac{v_p}{4+E^2}\left(
 * 4\cos\omega - 2E\sin\omega \right)
 * \f]
 * (where once again \f$4\cos\omega\f$ and \f$2\sin\omega\f$ can be precomputed).
 *
 * This routine does not account for relativistic timing variations, and
 * issues warnings or errors based on the criterea of
 * \eqref{eq_relativistic-orbit} in GenerateEllipticSpinOrbitCW().
 * The routine will also warn if
 * it seems likely that \c REAL8 precision may not be sufficient to
 * track the orbit accurately.  We estimate that numerical errors could
 * cause the number of computed wave cycles to vary by
 * \f[
 * \Delta N \lesssim f_0 T\epsilon\left[
 * \sim6+\ln\left(|C|+\sqrt{|C|^2+1}\right)\right] \;,
 * \f]
 * where \f$|C|\f$ is the maximum magnitude of the variable \f$C\f$ over the
 * course of the computation, \f$f_0T\f$ is the approximate total number of
 * wave cycles over the computation, and \f$\epsilon\approx2\times10^{-16}\f$
 * is the fractional precision of \c REAL8 arithmetic.  If this
 * estimate exceeds 0.01 cycles, a warning is issued.
 */
void
LALGenerateParabolicSpinOrbitCW( LALStatus             *stat,
				 PulsarCoherentGW            *output,
				 SpinOrbitCWParamStruc *params )
{
  UINT4 n, i;              /* number of and index over samples */
  UINT4 nSpin = 0, j;      /* number of and index over spindown terms */
  REAL8 t, dt, tPow;       /* time, interval, and t raised to a power */
  REAL8 phi0, f0, twopif0; /* initial phase, frequency, and 2*pi*f0 */
  REAL8 f, fPrev;      /* current and previous values of frequency */
  REAL4 df = 0.0;      /* maximum difference between f and fPrev */
  REAL8 phi;           /* current value of phase */
  REAL8 vp;            /* projected speed at periapsis */
  REAL8 argument;      /* argument of periapsis */
  REAL8 fourCosOmega;  /* four times the cosine of argument */
  REAL8 twoSinOmega;   /* two times the sine of argument */
  REAL8 vpCosOmega;    /* vp times cosine of argument */
  REAL8 vpSinOmega;    /* vp times sine of argument */
  REAL8 vpSinOmega2;   /* vpSinOmega squared */
  REAL8 vDot6;         /* 6 times angular speed at periapsis */
  REAL8 oneBy12vDot;   /* one over (12 times angular speed) */
  REAL8 pBy3;          /* constant sqrt(p/3) in cubic equation */
  REAL8 p32;           /* constant (p/3)^1.5 in cubic equation */
  REAL8 c, c0, dc;     /* C variable, offset, and step increment */
  REAL8 e, e2, e0;     /* E variable, E^2, and constant piece of E */
  REAL8 de;            /* coefficient of sinh() piece of E */
  REAL8 tpOff;         /* orbit epoch - time series epoch (s) */
  REAL8 spinOff;       /* spin epoch - orbit epoch (s) */
  REAL8 *fSpin = NULL; /* pointer to Taylor coefficients */
  REAL4 *fData;        /* pointer to frequency data */
  REAL8 *phiData;      /* pointer to phase data */

  INITSTATUS(stat);
  ATTATCHSTATUSPTR( stat );

  /* Make sure parameter and output structures exist. */
  ASSERT( params, stat, GENERATESPINORBITCWH_ENUL,
	  GENERATESPINORBITCWH_MSGENUL );
  ASSERT( output, stat, GENERATESPINORBITCWH_ENUL,
	  GENERATESPINORBITCWH_MSGENUL );

  /* Make sure output fields don't exist. */
  ASSERT( !( output->a ), stat, GENERATESPINORBITCWH_EOUT,
	  GENERATESPINORBITCWH_MSGEOUT );
  ASSERT( !( output->f ), stat, GENERATESPINORBITCWH_EOUT,
	  GENERATESPINORBITCWH_MSGEOUT );
  ASSERT( !( output->phi ), stat, GENERATESPINORBITCWH_EOUT,
	  GENERATESPINORBITCWH_MSGEOUT );
  ASSERT( !( output->shift ), stat, GENERATESPINORBITCWH_EOUT,
	  GENERATESPINORBITCWH_MSGEOUT );

  /* If Taylor coeficients are specified, make sure they exist. */
  if ( params->f ) {
    ASSERT( params->f->data, stat, GENERATESPINORBITCWH_ENUL,
	    GENERATESPINORBITCWH_MSGENUL );
    nSpin = params->f->length;
    fSpin = params->f->data;
  }

  /* Set up some constants (to avoid repeated calculation or
     dereferencing), and make sure they have acceptable values. */
  vp = params->rPeriNorm*params->angularSpeed;
  vDot6 = 6.0*params->angularSpeed;
  n = params->length;
  dt = params->deltaT;
  f0 = fPrev = params->f0;
  if ( params->oneMinusEcc != 0.0 ) {
    ABORT( stat, GENERATESPINORBITCWH_EECC,
	   GENERATESPINORBITCWH_MSGEECC );
  }
  if ( vp >= 1.0 ) {
    ABORT( stat, GENERATESPINORBITCWH_EFTL,
	   GENERATESPINORBITCWH_MSGEFTL );
  }
  if ( vp <= 0.0 || dt <= 0.0 || f0 <= 0.0 || vDot6 <= 0.0 ||
       n == 0 ) {
    ABORT( stat, GENERATESPINORBITCWH_ESGN,
	   GENERATESPINORBITCWH_MSGESGN );
  }
#ifndef NDEBUG
  if ( lalDebugLevel & LALWARNING ) {
    if ( f0*n*dt*vp*vp > 0.5 )
      LALWarning( stat, "Orbit may have significant relativistic"
		  " effects that are not included" );
  }
#endif

  /* Compute offset between time series epoch and periapsis, and
     betweem periapsis and spindown reference epoch. */
  tpOff = (REAL8)( params->orbitEpoch.gpsSeconds -
		   params->epoch.gpsSeconds );
  tpOff += 1.0e-9 * (REAL8)( params->orbitEpoch.gpsNanoSeconds -
			     params->epoch.gpsNanoSeconds );
  spinOff = (REAL8)( params->orbitEpoch.gpsSeconds -
		     params->spinEpoch.gpsSeconds );
  spinOff += 1.0e-9 * (REAL8)( params->orbitEpoch.gpsNanoSeconds -
			       params->spinEpoch.gpsNanoSeconds );

  /* Set up some other constants. */
  twopif0 = f0*LAL_TWOPI;
  phi0 = params->phi0;
  argument = params->omega;
  oneBy12vDot = 0.5/vDot6;
  fourCosOmega = 4.0*cos( argument );
  twoSinOmega = 2.0*sin( argument );
  vpCosOmega = 0.25*vp*fourCosOmega;
  vpSinOmega = 0.5*vp*twoSinOmega;
  vpSinOmega2 = vpSinOmega*vpSinOmega;
  pBy3 = sqrt( 4.0*( 1.0 + vpCosOmega ) - vpSinOmega2 );
  p32 = 1.0/( pBy3*pBy3*pBy3 );
  c0 = p32*( vpSinOmega*( 6.0*vpCosOmega - 12.0 + vpSinOmega2 ) -
	     tpOff*vDot6 );
  dc = p32*vDot6*dt;
  e0 = 3.0*vpSinOmega;
  de = 2.0*pBy3;

  /* Check whether REAL8 precision is good enough. */
#ifndef NDEBUG
  if ( lalDebugLevel & LALWARNING ) {
    REAL8 x = fabs( c0 + n*dc ); /* a temporary computation variable */
    if ( x < fabs( c0 ) )
      x = fabs( c0 );
    x = 6.0 + log( x + sqrt( x*x + 1.0 ) );
    if ( LAL_REAL8_EPS*f0*dt*n*x > 0.01 )
      LALWarning( stat, "REAL8 arithmetic may not have sufficient"
		  " precision for this orbit" );
  }
#endif

  /* Allocate output structures. */
  if ( ( output->a = (REAL4TimeVectorSeries *)
	 LALMalloc( sizeof(REAL4TimeVectorSeries) ) ) == NULL ) {
    ABORT( stat, GENERATESPINORBITCWH_EMEM,
	   GENERATESPINORBITCWH_MSGEMEM );
  }
  memset( output->a, 0, sizeof(REAL4TimeVectorSeries) );
  if ( ( output->f = (REAL4TimeSeries *)
	 LALMalloc( sizeof(REAL4TimeSeries) ) ) == NULL ) {
    LALFree( output->a ); output->a = NULL;
    ABORT( stat, GENERATESPINORBITCWH_EMEM,
	   GENERATESPINORBITCWH_MSGEMEM );
  }
  memset( output->f, 0, sizeof(REAL4TimeSeries) );
  if ( ( output->phi = (REAL8TimeSeries *)
	 LALMalloc( sizeof(REAL8TimeSeries) ) ) == NULL ) {
    LALFree( output->a ); output->a = NULL;
    LALFree( output->f ); output->f = NULL;
    ABORT( stat, GENERATESPINORBITCWH_EMEM,
	   GENERATESPINORBITCWH_MSGEMEM );
  }
  memset( output->phi, 0, sizeof(REAL8TimeSeries) );

  /* Set output structure metadata fields. */
  output->position = params->position;
  output->psi = params->psi;
  output->a->epoch = output->f->epoch = output->phi->epoch
    = params->epoch;
  output->a->deltaT = n*params->deltaT;
  output->f->deltaT = output->phi->deltaT = params->deltaT;
  output->a->sampleUnits = lalStrainUnit;
  output->f->sampleUnits = lalHertzUnit;
  output->phi->sampleUnits = lalDimensionlessUnit;
  snprintf( output->a->name, LALNameLength, "CW amplitudes" );
  snprintf( output->f->name, LALNameLength, "CW frequency" );
  snprintf( output->phi->name, LALNameLength, "CW phase" );

  /* Allocate phase and frequency arrays. */
  LALSCreateVector( stat->statusPtr, &( output->f->data ), n );
  BEGINFAIL( stat ) {
    LALFree( output->a );   output->a = NULL;
    LALFree( output->f );   output->f = NULL;
    LALFree( output->phi ); output->phi = NULL;
  } ENDFAIL( stat );
  LALDCreateVector( stat->statusPtr, &( output->phi->data ), n );
  BEGINFAIL( stat ) {
    TRY( LALSDestroyVector( stat->statusPtr, &( output->f->data ) ),
	 stat );
    LALFree( output->a );   output->a = NULL;
    LALFree( output->f );   output->f = NULL;
    LALFree( output->phi ); output->phi = NULL;
  } ENDFAIL( stat );

  /* Allocate and fill amplitude array. */
  {
    CreateVectorSequenceIn in; /* input to create output->a */
    in.length = 2;
    in.vectorLength = 2;
    LALSCreateVectorSequence( stat->statusPtr, &(output->a->data), &in );
    BEGINFAIL( stat ) {
      TRY( LALSDestroyVector( stat->statusPtr, &( output->f->data ) ),
	   stat );
      TRY( LALDDestroyVector( stat->statusPtr, &( output->phi->data ) ),
	   stat );
      LALFree( output->a );   output->a = NULL;
      LALFree( output->f );   output->f = NULL;
      LALFree( output->phi ); output->phi = NULL;
    } ENDFAIL( stat );
    output->a->data->data[0] = output->a->data->data[2] = params->aPlus;
    output->a->data->data[1] = output->a->data->data[3] = params->aCross;
  }

  /* Fill frequency and phase arrays. */
  fData = output->f->data->data;
  phiData = output->phi->data->data;
  for ( i = 0; i < n; i++ ) {

    /* Compute emission time. */
    c = c0 + dc*i;
    if ( c > 0 )
      e = e0 + de*sinh( log( c + sqrt( c*c + 1.0 ) )/3.0 );
    else
      e = e0 + de*sinh( -log( -c + sqrt( c*c + 1.0 ) )/3.0 );
    e2 = e*e;
    phi = t = tPow = oneBy12vDot*e*( 12.0 + e2 );

    /* Compute source emission phase and frequency. */
    f = 1.0;
    for ( j = 0; j < nSpin; j++ ) {
      f += fSpin[j]*tPow;
      phi += fSpin[j]*( tPow*=t )/( j + 2.0 );
    }

    /* Appy frequency Doppler shift. */
    f *= f0 / ( 1.0 + vp*( fourCosOmega - e*twoSinOmega )
		/( 4.0 + e2 ) );
    phi *= twopif0;
    if ( fabs( f - fPrev ) > df )
      df = fabs( f - fPrev );
    *(fData++) = fPrev = f;
    *(phiData++) = phi + phi0;
  }

  /* Set output field and return. */
  params->dfdt = df*dt;
  DETATCHSTATUSPTR( stat );
  RETURN( stat );
}
Exemplo n.º 2
0
void
LALGenerateRing(
    LALStatus          *stat,
    CoherentGW         *output,
    REAL4TimeSeries    UNUSED *series,
    SimRingdownTable   *simRingdown,
    RingParamStruc     *params
    )

{
  UINT4 i;      /* number of and index over samples */
  REAL8 t, dt;         /* time, interval */
  REAL8 gtime ;    /* central time, decay time */
  REAL8 f0, quality;   /* frequency and quality factor */
  REAL8 twopif0;       /* 2*pi*f0 */
  REAL4 h0;            /* peak strain for ringdown */
  REAL4 *fData;        /* pointer to frequency data */
  REAL8 *phiData;      /* pointer to phase data */
  REAL8 init_phase;    /*initial phase of injection */
  REAL4 *aData;        /* pointer to frequency data */
  UINT4 nPointInj; /* number of data points in a block */
#if 0
  UINT4 n;
  REAL8 t0;
  REAL4TimeSeries signalvec; /* start time of block that injection is injected into */
  LALTimeInterval  interval;
  INT8 geoc_tns;       /* geocentric_start_time of the injection in ns */
  INT8 block_tns;      /* start time of block in ns */
  REAL8 deltaTns;
  INT8 inj_diff;       /* time between start of segment and injection */
  LALTimeInterval dummyInterval;
#endif

  INITSTATUS(stat);
  ATTATCHSTATUSPTR( stat );

  /* Make sure parameter and output structures exist. */
  ASSERT( params, stat, GENERATERINGH_ENUL,
	  GENERATERINGH_MSGENUL );
  ASSERT( output, stat, GENERATERINGH_ENUL,
	  GENERATERINGH_MSGENUL );

  /* Make sure output fields don't exist. */
  ASSERT( !( output->a ), stat, GENERATERINGH_EOUT,
	  GENERATERINGH_MSGEOUT );
  ASSERT( !( output->f ), stat, GENERATERINGH_EOUT,
	  GENERATERINGH_MSGEOUT );
  ASSERT( !( output->phi ), stat, GENERATERINGH_EOUT,
	  GENERATERINGH_MSGEOUT );
  ASSERT( !( output->shift ), stat, GENERATERINGH_EOUT,
	  GENERATERINGH_MSGEOUT );


  /* Set up some other constants, to avoid repeated dereferencing. */
  dt = params->deltaT;
/* N_point = 2 * floor(0.5+ 1/ dt); */

  nPointInj = 163840;

  /* Generic ring parameters */
  h0 = simRingdown->amplitude;
  quality = (REAL8)simRingdown->quality;
  f0 = (REAL8)simRingdown->frequency;
  twopif0 = f0*LAL_TWOPI;
  init_phase = simRingdown->phase;


  if ( ( output->a = (REAL4TimeVectorSeries *)
	 LALMalloc( sizeof(REAL4TimeVectorSeries) ) ) == NULL ) {
    ABORT( stat, GENERATERINGH_EMEM, GENERATERINGH_MSGEMEM );
  }
  memset( output->a, 0, sizeof(REAL4TimeVectorSeries) );
  if ( ( output->f = (REAL4TimeSeries *)
	 LALMalloc( sizeof(REAL4TimeSeries) ) ) == NULL ) {
    LALFree( output->a ); output->a = NULL;
    ABORT( stat, GENERATERINGH_EMEM, GENERATERINGH_MSGEMEM );
  }
  memset( output->f, 0, sizeof(REAL4TimeSeries) );
  if ( ( output->phi = (REAL8TimeSeries *)
	 LALMalloc( sizeof(REAL8TimeSeries) ) ) == NULL ) {
    LALFree( output->a ); output->a = NULL;
    LALFree( output->f ); output->f = NULL;
    ABORT( stat, GENERATERINGH_EMEM, GENERATERINGH_MSGEMEM );
  }
  memset( output->phi, 0, sizeof(REAL8TimeSeries) );

  /* Set output structure metadata fields. */
  output->position.longitude = simRingdown->longitude;
  output->position.latitude = simRingdown->latitude;
  output->position.system = params->system;
  output->psi = simRingdown->polarization;
   /* set epoch of output time series to that of the block */
  output->a->epoch = output->f->epoch = output->phi->epoch = simRingdown->geocent_start_time;
  output->a->deltaT = params->deltaT;
  output->f->deltaT = output->phi->deltaT = params->deltaT;
  output->a->sampleUnits = lalStrainUnit;
  output->f->sampleUnits = lalHertzUnit;
  output->phi->sampleUnits = lalDimensionlessUnit;
  snprintf( output->a->name, LALNameLength, "Ring amplitudes" );
  snprintf( output->f->name, LALNameLength, "Ring frequency" );
  snprintf( output->phi->name, LALNameLength, "Ring phase" );


  /* Allocate phase and frequency arrays. */
  LALSCreateVector( stat->statusPtr, &( output->f->data ), nPointInj );
  BEGINFAIL( stat ) {
    LALFree( output->a );   output->a = NULL;
    LALFree( output->f );   output->f = NULL;
    LALFree( output->phi ); output->phi = NULL;
  } ENDFAIL( stat );

  LALDCreateVector( stat->statusPtr, &( output->phi->data ), nPointInj );
  BEGINFAIL( stat ) {
    TRY( LALSDestroyVector( stat->statusPtr, &( output->f->data ) ),
	 stat );
    LALFree( output->a );   output->a = NULL;
    LALFree( output->f );   output->f = NULL;
    LALFree( output->phi ); output->phi = NULL;
  } ENDFAIL( stat );


  /* Allocate amplitude array. */
  {
    CreateVectorSequenceIn in; /* input to create output->a */
    in.length = nPointInj;
    in.vectorLength = 2;
    LALSCreateVectorSequence( stat->statusPtr, &(output->a->data), &in );
    BEGINFAIL( stat ) {
      TRY( LALSDestroyVector( stat->statusPtr, &( output->f->data ) ),
	   stat );
      TRY( LALDDestroyVector( stat->statusPtr, &( output->phi->data ) ),
	   stat );
      LALFree( output->a );   output->a = NULL;
      LALFree( output->f );   output->f = NULL;
      LALFree( output->phi ); output->phi = NULL;
    } ENDFAIL( stat );
  }


  /*  set arrays to zero */
  memset( output->f->data->data, 0, sizeof( REAL4 ) *  output->f->data->length );
  memset( output->phi->data->data, 0, sizeof( REAL8 ) * output->phi->data->length );
  memset( output->a->data->data, 0, sizeof( REAL4 ) *
      output->a->data->length * output->a->data->vectorLength );

/* Fill frequency and phase arrays starting at time of injection NOT start */
  fData = output->f->data->data;
  phiData = output->phi->data->data;
  aData = output->a->data->data;

  if ( !( strcmp( simRingdown->waveform, "Ringdown" ) ) )
  {
    for ( i = 0; i < nPointInj; i++ )
    {
      t = i * dt;
      gtime = twopif0 / 2 / quality * t ;
      *(fData++)   = f0;
      *(phiData++) = twopif0 * t+init_phase;
      *(aData++) = h0 * ( 1.0 + pow( cos( simRingdown->inclination ), 2 ) ) *
        exp( - gtime );
      *(aData++) = h0* 2.0 * cos( simRingdown->inclination ) * exp( - gtime );
    }
  }
  else
  {
    ABORT( stat, GENERATERINGH_ETYP, GENERATERINGH_MSGETYP );
  }


/* Set output field and return. */
  DETATCHSTATUSPTR( stat );
  RETURN( stat );
}
/**
 * \author Creighton, T. D.
 *
 * \brief Computes a continuous waveform with frequency drift and Doppler
 * modulation from a hyperbolic orbital trajectory.
 *
 * This function computes a quaiperiodic waveform using the spindown and
 * orbital parameters in <tt>*params</tt>, storing the result in
 * <tt>*output</tt>.
 *
 * In the <tt>*params</tt> structure, the routine uses all the "input"
 * fields specified in \ref GenerateSpinOrbitCW_h, and sets all of the
 * "output" fields.  If <tt>params-\>f</tt>=\c NULL, no spindown
 * modulation is performed.  If <tt>params-\>oneMinusEcc</tt>\f$\not<0\f$ (a
 * non-hyperbolic orbit), or if
 * <tt>params-\>rPeriNorm</tt>\f$\times\f$<tt>params-\>angularSpeed</tt>\f$\geq1\f$
 * (faster-than-light speed at periapsis), an error is returned.
 *
 * In the <tt>*output</tt> structure, the field <tt>output-\>h</tt> is
 * ignored, but all other pointer fields must be set to \c NULL.  The
 * function will create and allocate space for <tt>output-\>a</tt>,
 * <tt>output-\>f</tt>, and <tt>output-\>phi</tt> as necessary.  The
 * <tt>output-\>shift</tt> field will remain set to \c NULL.
 *
 * ### Algorithm ###
 *
 * For hyperbolic orbits, we combine \eqref{eq_spinorbit-tr},
 * \eqref{eq_spinorbit-t}, and \eqref{eq_spinorbit-upsilon} to get \f$t_r\f$
 * directly as a function of \f$E\f$:
 * \f{eqnarray}{
 * \label{eq_tr-e3}
 * t_r = t_p & + & \left(\frac{r_p \sin i}{c}\right)\sin\omega\\
 * & + & \frac{1}{n} \left( -E +
 * \left[v_p(e-1)\cos\omega + e\right]\sinh E
 * - \left[v_p\sqrt{\frac{e-1}{e+1}}\sin\omega\right]
 * [\cosh E - 1]\right) \;,
 * \f}
 * where \f$v_p=r_p\dot{\upsilon}_p\sin i/c\f$ is a normalized velocity at
 * periapsis and \f$n=\dot{\upsilon}_p\sqrt{(1-e)^3/(1+e)}\f$ is a normalized
 * angular speed for the orbit (the hyperbolic analogue of the mean
 * angular speed for closed orbits).  For simplicity we write this as:
 * \f{equation}{
 * \label{eq_tr-e4}
 * t_r = T_p + \frac{1}{n}\left( E + A\sinh E + B[\cosh E - 1] \right) \;,
 * \f}
 *
 * \figure{inject_hanomaly,eps,0.23,Function to be inverted to find eccentric anomaly}
 *
 * where \f$T_p\f$ is the \e observed time of periapsis passage.  Thus
 * the key numerical procedure in this routine is to invert the
 * expression \f$x=E+A\sinh E+B(\cosh E - 1)\f$ to get \f$E(x)\f$.  This function
 * is sketched to the right (solid line), along with an approximation
 * used for making an initial guess (dotted line), as described later.
 *
 * We note that \f$A^2-B^2<1\f$, although it approaches 1 when
 * \f$e\rightarrow1\f$, or when \f$v_p\rightarrow1\f$ and either \f$e=0\f$ or
 * \f$\omega=\pi\f$.  Except in this limit, Newton-Raphson methods will
 * converge rapidly for any initial guess.  In this limit, though, the
 * slope \f$dx/dE\f$ approaches zero at \f$E=0\f$, and an initial guess or
 * iteration landing near this point will send the next iteration off to
 * unacceptably large or small values.  A hybrid root-finding strategy is
 * used to deal with this, and with the exponential behaviour of \f$x\f$ at
 * large \f$E\f$.
 *
 * First, we compute \f$x=x_{\pm1}\f$ at \f$E=\pm1\f$.  If the desired \f$x\f$ lies
 * in this range, we use a straightforward Newton-Raphson root finder,
 * with the constraint that all guesses of \f$E\f$ are restricted to the
 * domain \f$[-1,1]\f$.  This guarantees that the scheme will eventually find
 * itself on a uniformly-convergent trajectory.
 *
 * Second, for \f$E\f$ outside of this range, \f$x\f$ is dominated by the
 * exponential terms: \f$x\approx\frac{1}{2}(A+B)\exp(E)\f$ for \f$E\gg1\f$, and
 * \f$x\approx-\frac{1}{2}(A-B)\exp(-E)\f$ for \f$E\ll-1\f$.  We therefore do an
 * \e approximate Newton-Raphson iteration on the function \f$\ln|x|\f$,
 * where the approximation is that we take \f$d\ln|x|/d|E|\approx1\f$.  This
 * involves computing an extra logarithm inside the loop, but gives very
 * rapid convergence to high precision, since \f$\ln|x|\f$ is very nearly
 * linear in these regions.
 *
 * At the start of the algorithm, we use an initial guess of
 * \f$E=-\ln[-2(x-x_{-1})/(A-B)-\exp(1)]\f$ for \f$x<x_{-1}\f$, \f$E=x/x_{-1}\f$ for
 * \f$x_{-1}\leq x\leq0\f$, \f$E=x/x_{+1}\f$ for \f$0\leq x\leq x_{+1}\f$, or
 * \f$E=\ln[2(x-x_{+1})/(A+B)-\exp(1)]\f$ for \f$x>x_{+1}\f$.  We refine this
 * guess until we get agreement to within 0.01 parts in part in
 * \f$N_\mathrm{cyc}\f$ (where \f$N_\mathrm{cyc}\f$ is the larger of the number
 * of wave cycles in a time \f$2\pi/n\f$, or the number of wave cycles in the
 * entire waveform being generated), or one part in \f$10^{15}\f$ (an order
 * of magnitude off the best precision possible with \c REAL8
 * numbers).  The latter case indicates that \c REAL8 precision may
 * fail to give accurate phasing, and one should consider modeling the
 * orbit as a set of Taylor frequency coefficients \'{a} la
 * <tt>LALGenerateTaylorCW()</tt>.  On subsequent timesteps, we use the
 * previous timestep as an initial guess, which is good so long as the
 * timesteps are much smaller than \f$1/n\f$.
 *
 * Once a value of \f$E\f$ is found for a given timestep in the output
 * series, we compute the system time \f$t\f$ via \eqref{eq_spinorbit-t},
 * and use it to determine the wave phase and (non-Doppler-shifted)
 * frequency via \eqref{eq_taylorcw-freq}
 * and \eqref{eq_taylorcw-phi}.  The Doppler shift on the frequency is
 * then computed using \eqref{eq_spinorbit-upsilon}
 * and \eqref{eq_orbit-rdot}.  We use \f$\upsilon\f$ as an intermediate in
 * the Doppler shift calculations, since expressing \f$\dot{R}\f$ directly in
 * terms of \f$E\f$ results in expression of the form \f$(e-1)/(e\cosh E-1)\f$,
 * which are difficult to simplify and face precision losses when
 * \f$E\sim0\f$ and \f$e\rightarrow1\f$.  By contrast, solving for \f$\upsilon\f$ is
 * numerically stable provided that the system <tt>atan2()</tt> function is
 * well-designed.
 *
 * This routine does not account for relativistic timing variations, and
 * issues warnings or errors based on the criterea of
 * \eqref{eq_relativistic-orbit} in \ref LALGenerateEllipticSpinOrbitCW().
 */
void
LALGenerateHyperbolicSpinOrbitCW( LALStatus             *stat,
				  PulsarCoherentGW            *output,
				  SpinOrbitCWParamStruc *params )
{
  UINT4 n, i;              /* number of and index over samples */
  UINT4 nSpin = 0, j;      /* number of and index over spindown terms */
  REAL8 t, dt, tPow;       /* time, interval, and t raised to a power */
  REAL8 phi0, f0, twopif0; /* initial phase, frequency, and 2*pi*f0 */
  REAL8 f, fPrev;          /* current and previous values of frequency */
  REAL4 df = 0.0;          /* maximum difference between f and fPrev */
  REAL8 phi;               /* current value of phase */
  REAL8 vDotAvg;           /* nomalized orbital angular speed */
  REAL8 vp;                /* projected speed at periapsis */
  REAL8 upsilon, argument; /* true anomaly, and argument of periapsis */
  REAL8 eCosOmega;         /* eccentricity * cosine of argument */
  REAL8 tPeriObs;          /* time of observed periapsis */
  REAL8 spinOff;           /* spin epoch - orbit epoch */
  REAL8 x;                 /* observed ``mean anomaly'' */
  REAL8 xPlus, xMinus;     /* limits where exponentials dominate */
  REAL8 dx, dxMax;         /* current and target errors in x */
  REAL8 a, b;              /* constants in equation for x */
  REAL8 ecc;               /* orbital eccentricity */
  REAL8 eccMinusOne, eccPlusOne; /* ecc - 1 and ecc + 1 */
  REAL8 e;                       /* hyperbolic anomaly */
  REAL8 sinhe, coshe;            /* sinh of e, and cosh of e minus 1 */
  REAL8 *fSpin = NULL;           /* pointer to Taylor coefficients */
  REAL4 *fData;                  /* pointer to frequency data */
  REAL8 *phiData;                /* pointer to phase data */

  INITSTATUS(stat);
  ATTATCHSTATUSPTR( stat );

  /* Make sure parameter and output structures exist. */
  ASSERT( params, stat, GENERATESPINORBITCWH_ENUL,
	  GENERATESPINORBITCWH_MSGENUL );
  ASSERT( output, stat, GENERATESPINORBITCWH_ENUL,
	  GENERATESPINORBITCWH_MSGENUL );

  /* Make sure output fields don't exist. */
  ASSERT( !( output->a ), stat, GENERATESPINORBITCWH_EOUT,
	  GENERATESPINORBITCWH_MSGEOUT );
  ASSERT( !( output->f ), stat, GENERATESPINORBITCWH_EOUT,
	  GENERATESPINORBITCWH_MSGEOUT );
  ASSERT( !( output->phi ), stat, GENERATESPINORBITCWH_EOUT,
	  GENERATESPINORBITCWH_MSGEOUT );
  ASSERT( !( output->shift ), stat, GENERATESPINORBITCWH_EOUT,
	  GENERATESPINORBITCWH_MSGEOUT );

  /* If Taylor coeficients are specified, make sure they exist. */
  if ( params->f ) {
    ASSERT( params->f->data, stat, GENERATESPINORBITCWH_ENUL,
	    GENERATESPINORBITCWH_MSGENUL );
    nSpin = params->f->length;
    fSpin = params->f->data;
  }

  /* Set up some constants (to avoid repeated calculation or
     dereferencing), and make sure they have acceptable values. */
  eccMinusOne = -params->oneMinusEcc;
  ecc = 1.0 + eccMinusOne;
  eccPlusOne = 2.0 + eccMinusOne;
  if ( eccMinusOne <= 0.0 ) {
    ABORT( stat, GENERATESPINORBITCWH_EECC,
	   GENERATESPINORBITCWH_MSGEECC );
  }
  vp = params->rPeriNorm*params->angularSpeed;
  vDotAvg = params->angularSpeed
    *sqrt( eccMinusOne*eccMinusOne*eccMinusOne/eccPlusOne );
  n = params->length;
  dt = params->deltaT;
  f0 = fPrev = params->f0;
  if ( vp >= 1.0 ) {
    ABORT( stat, GENERATESPINORBITCWH_EFTL,
	   GENERATESPINORBITCWH_MSGEFTL );
  }
  if ( vp <= 0.0 || dt <= 0.0 || f0 <= 0.0 || vDotAvg <= 0.0 ||
       n == 0 ) {
    ABORT( stat, GENERATESPINORBITCWH_ESGN,
	   GENERATESPINORBITCWH_MSGESGN );
  }

  /* Set up some other constants. */
  twopif0 = f0*LAL_TWOPI;
  phi0 = params->phi0;
  argument = params->omega;
  a = vp*eccMinusOne*cos( argument ) + ecc;
  b = -vp*sqrt( eccMinusOne/eccPlusOne )*sin( argument );
  eCosOmega = ecc*cos( argument );
  if ( n*dt*vDotAvg > LAL_TWOPI )
    dxMax = 0.01/( f0*n*dt );
  else
    dxMax = 0.01/( f0*LAL_TWOPI/vDotAvg );
  if ( dxMax < 1.0e-15 ) {
    dxMax = 1.0e-15;
#ifndef NDEBUG
    LALWarning( stat, "REAL8 arithmetic may not have sufficient"
		" precision for this orbit" );
#endif
  }
#ifndef NDEBUG
  if ( lalDebugLevel & LALWARNING ) {
    REAL8 tau = n*dt;
    if ( tau > LAL_TWOPI/vDotAvg )
      tau = LAL_TWOPI/vDotAvg;
    if ( f0*tau*vp*vp*ecc/eccPlusOne > 0.25 )
      LALWarning( stat, "Orbit may have significant relativistic"
		  " effects that are not included" );
  }
#endif

  /* Compute offset between time series epoch and observed periapsis,
     and betweem true periapsis and spindown reference epoch. */
  tPeriObs = (REAL8)( params->orbitEpoch.gpsSeconds -
		      params->epoch.gpsSeconds );
  tPeriObs += 1.0e-9 * (REAL8)( params->orbitEpoch.gpsNanoSeconds -
				params->epoch.gpsNanoSeconds );
  tPeriObs += params->rPeriNorm*sin( params->omega );
  spinOff = (REAL8)( params->orbitEpoch.gpsSeconds -
		     params->spinEpoch.gpsSeconds );
  spinOff += 1.0e-9 * (REAL8)( params->orbitEpoch.gpsNanoSeconds -
			       params->spinEpoch.gpsNanoSeconds );

  /* Determine bounds of hybrid root-finding algorithm, and initial
     guess for e. */
  xMinus = 1.0 + a*sinh( -1.0 ) + b*cosh( -1.0 ) - b;
  xPlus = -1.0 + a*sinh( 1.0 ) + b*cosh( 1.0 ) - b;
  x = -vDotAvg*tPeriObs;
  if ( x < xMinus )
    e = -log( -2.0*( x - xMinus )/( a - b ) - exp( 1.0 ) );
  else if ( x <= 0 )
    e = x/xMinus;
  else if ( x <= xPlus )
    e = x/xPlus;
  else
    e = log( 2.0*( x - xPlus )/( a + b ) - exp( 1.0 ) );
  sinhe = sinh( e );
  coshe = cosh( e ) - 1.0;

  /* Allocate output structures. */
  if ( ( output->a = (REAL4TimeVectorSeries *)
	 LALMalloc( sizeof(REAL4TimeVectorSeries) ) ) == NULL ) {
    ABORT( stat, GENERATESPINORBITCWH_EMEM,
	   GENERATESPINORBITCWH_MSGEMEM );
  }
  memset( output->a, 0, sizeof(REAL4TimeVectorSeries) );
  if ( ( output->f = (REAL4TimeSeries *)
	 LALMalloc( sizeof(REAL4TimeSeries) ) ) == NULL ) {
    LALFree( output->a ); output->a = NULL;
    ABORT( stat, GENERATESPINORBITCWH_EMEM,
	   GENERATESPINORBITCWH_MSGEMEM );
  }
  memset( output->f, 0, sizeof(REAL4TimeSeries) );
  if ( ( output->phi = (REAL8TimeSeries *)
	 LALMalloc( sizeof(REAL8TimeSeries) ) ) == NULL ) {
    LALFree( output->a ); output->a = NULL;
    LALFree( output->f ); output->f = NULL;
    ABORT( stat, GENERATESPINORBITCWH_EMEM,
	   GENERATESPINORBITCWH_MSGEMEM );
  }
  memset( output->phi, 0, sizeof(REAL8TimeSeries) );

  /* Set output structure metadata fields. */
  output->position = params->position;
  output->psi = params->psi;
  output->a->epoch = output->f->epoch = output->phi->epoch
    = params->epoch;
  output->a->deltaT = n*params->deltaT;
  output->f->deltaT = output->phi->deltaT = params->deltaT;
  output->a->sampleUnits = lalStrainUnit;
  output->f->sampleUnits = lalHertzUnit;
  output->phi->sampleUnits = lalDimensionlessUnit;
  snprintf( output->a->name, LALNameLength, "CW amplitudes" );
  snprintf( output->f->name, LALNameLength, "CW frequency" );
  snprintf( output->phi->name, LALNameLength, "CW phase" );

  /* Allocate phase and frequency arrays. */
  LALSCreateVector( stat->statusPtr, &( output->f->data ), n );
  BEGINFAIL( stat ) {
    LALFree( output->a );   output->a = NULL;
    LALFree( output->f );   output->f = NULL;
    LALFree( output->phi ); output->phi = NULL;
  } ENDFAIL( stat );
  LALDCreateVector( stat->statusPtr, &( output->phi->data ), n );
  BEGINFAIL( stat ) {
    TRY( LALSDestroyVector( stat->statusPtr, &( output->f->data ) ),
	 stat );
    LALFree( output->a );   output->a = NULL;
    LALFree( output->f );   output->f = NULL;
    LALFree( output->phi ); output->phi = NULL;
  } ENDFAIL( stat );

  /* Allocate and fill amplitude array. */
  {
    CreateVectorSequenceIn in; /* input to create output->a */
    in.length = 2;
    in.vectorLength = 2;
    LALSCreateVectorSequence( stat->statusPtr, &(output->a->data), &in );
    BEGINFAIL( stat ) {
      TRY( LALSDestroyVector( stat->statusPtr, &( output->f->data ) ),
	   stat );
      TRY( LALDDestroyVector( stat->statusPtr, &( output->phi->data ) ),
	   stat );
      LALFree( output->a );   output->a = NULL;
      LALFree( output->f );   output->f = NULL;
      LALFree( output->phi ); output->phi = NULL;
    } ENDFAIL( stat );
    output->a->data->data[0] = output->a->data->data[2] = params->aPlus;
    output->a->data->data[1] = output->a->data->data[3] = params->aCross;
  }

  /* Fill frequency and phase arrays. */
  fData = output->f->data->data;
  phiData = output->phi->data->data;
  for ( i = 0; i < n; i++ ) {

    x = vDotAvg*( i*dt - tPeriObs );

    /* Use approximate Newton-Raphson method on ln|x| if |x| > 1. */
    if ( x < xMinus ) {
      x = log( -x );
      while ( fabs( dx = log( e - a*sinhe - b*coshe ) - x ) > dxMax ) {
	e += dx;
	sinhe = sinh( e );
	coshe = cosh( e ) - 1.0;
      }
    }
    else if ( x > xPlus ) {
      x = log( x );
      while ( fabs( dx = log( -e + a*sinhe + b*coshe ) - x ) > dxMax ) {
	e -= dx;
	sinhe = sinh( e );
	coshe = cosh( e ) - 1.0;
      }
    }

    /* Use ordinary Newton-Raphson method on x if |x| <= 1. */
    else {
      while ( fabs( dx = -e + a*sinhe + b*coshe - x ) > dxMax ) {
	e -= dx/( -1.0 + a*coshe + a + b*sinhe );
	if ( e < -1.0 )
	  e = -1.0;
	else if ( e > 1.0 )
	  e = 1.0;
	sinhe = sinh( e );
	coshe = cosh( e ) - 1.0;
      }
    }

    /* Compute source emission time, phase, and frequency. */
    phi = t = tPow =
      ( ecc*sinhe - e )/vDotAvg + spinOff;
    f = 1.0;
    for ( j = 0; j < nSpin; j++ ) {
      f += fSpin[j]*tPow;
      phi += fSpin[j]*( tPow*=t )/( j + 2.0 );
    }

    /* Appy frequency Doppler shift. */
    upsilon = 2.0*atan2( sqrt( eccPlusOne*coshe ),
			 sqrt( eccMinusOne*( coshe + 2.0 ) ) );
    f *= f0 / ( 1.0 + vp*( cos( argument + upsilon ) + eCosOmega )
		/eccPlusOne );
    phi *= twopif0;
    if ( fabs( f - fPrev ) > df )
      df = fabs( f - fPrev );
    *(fData++) = fPrev = f;
    *(phiData++) = phi + phi0;
  }

  /* Set output field and return. */
  params->dfdt = df*dt;
  DETATCHSTATUSPTR( stat );
  RETURN( stat );
}
void
LALInspiralEccentricityForInjection(
			     LALStatus        *status,
			     CoherentGW       *waveform,
			     InspiralTemplate *params,
			     PPNParamStruc  *ppnParams
			     )
{

  INT4        count, i;
  REAL8       p, phiC;

  REAL4Vector a;           /* pointers to generated amplitude  data */
  REAL4Vector ff;          /* pointers to generated  frequency data */
  REAL8Vector phi;         /* generated phase data */

  CreateVectorSequenceIn in;

  CHAR message[256];

  InspiralInit paramsInit;

  INITSTATUS(status);
  ATTATCHSTATUSPTR(status);

  /* Make sure parameter and waveform structures exist. */
  ASSERT( params, status, LALINSPIRALH_ENULL, LALINSPIRALH_MSGENULL );
  ASSERT(waveform, status, LALINSPIRALH_ENULL, LALINSPIRALH_MSGENULL);
  ASSERT( !( waveform->a ), status, LALINSPIRALH_ENULL, LALINSPIRALH_MSGENULL );
  ASSERT( !( waveform->f ), status, LALINSPIRALH_ENULL, LALINSPIRALH_MSGENULL );
  ASSERT( !( waveform->phi ), status, LALINSPIRALH_ENULL, LALINSPIRALH_MSGENULL );

  /* Compute some parameters*/
  LALInspiralInit(status->statusPtr, params, &paramsInit);
  CHECKSTATUSPTR(status);

  if (paramsInit.nbins == 0){
      DETATCHSTATUSPTR(status);
      RETURN (status);
  }

  /* Now we can allocate memory and vector for coherentGW structure*/

  ff.length  = paramsInit.nbins;
  a.length   = 2* paramsInit.nbins;
  phi.length = paramsInit.nbins;

  ff.data = (REAL4 *) LALCalloc(paramsInit.nbins, sizeof(REAL4));
  a.data  = (REAL4 *) LALCalloc(2 * paramsInit.nbins, sizeof(REAL4));
  phi.data= (REAL8 *) LALCalloc(paramsInit.nbins, sizeof(REAL8));

  /* Check momory allocation is okay */
  if (!(ff.data) || !(a.data) || !(phi.data))
  {
    if (ff.data)  LALFree(ff.data);
    if (a.data)   LALFree(a.data);
    if (phi.data) LALFree(phi.data);

    ABORT( status, LALINSPIRALH_EMEM, LALINSPIRALH_MSGEMEM );
  }

  count = 0;

  /* Call the engine function */
  LALInspiralEccentricityEngine(status->statusPtr, NULL, NULL, &a, &ff, &phi, &count, params);
  BEGINFAIL( status )
  {
    LALFree(ff.data);
    LALFree(a.data);
    LALFree(phi.data);
  }
  ENDFAIL( status );

  p = phi.data[count-1];

  params->fFinal = ff.data[count-1];
  sprintf(message, "cycles = %f", p/(double)LAL_TWOPI);
  LALInfo(status, message);

  if ( (INT4)(p/LAL_TWOPI) < 2 ){
    sprintf(message, "The waveform has only %f cycles; we don't keep waveform with less than 2 cycles.",
	       p/(double)LAL_TWOPI );
    XLALPrintError("%s", message);
    LALWarning(status, message);
  }

      /*wrap the phase vector*/
      phiC =  phi.data[count-1] ;
      for (i = 0; i < count; i++)
	{
	  phi.data[i] =  phi.data[i] - phiC + ppnParams->phi;
	}

      /* Allocate the waveform structures. */
      if ( ( waveform->a = (REAL4TimeVectorSeries *)
	     LALCalloc(1, sizeof(REAL4TimeVectorSeries) ) ) == NULL ) {
	ABORT( status, LALINSPIRALH_EMEM,
	       LALINSPIRALH_MSGEMEM );
      }
      if ( ( waveform->f = (REAL4TimeSeries *)
	     LALCalloc(1, sizeof(REAL4TimeSeries) ) ) == NULL ) {
	LALFree( waveform->a ); waveform->a = NULL;
	ABORT( status, LALINSPIRALH_EMEM,
	       LALINSPIRALH_MSGEMEM );
      }
      if ( ( waveform->phi = (REAL8TimeSeries *)
	     LALCalloc(1, sizeof(REAL8TimeSeries) ) ) == NULL ) {
	LALFree( waveform->a ); waveform->a = NULL;
	LALFree( waveform->f ); waveform->f = NULL;
	ABORT( status, LALINSPIRALH_EMEM,
	       LALINSPIRALH_MSGEMEM );
      }


      in.length = (UINT4)(count);
      in.vectorLength = 2;

      LALSCreateVectorSequence( status->statusPtr, &( waveform->a->data ), &in );
      CHECKSTATUSPTR(status);

      LALSCreateVector( status->statusPtr, &( waveform->f->data ), count);
      CHECKSTATUSPTR(status);

      LALDCreateVector( status->statusPtr, &( waveform->phi->data ), count );
      CHECKSTATUSPTR(status);

      memcpy(waveform->f->data->data , ff.data, count*(sizeof(REAL4)));
      memcpy(waveform->a->data->data , a.data, 2*count*(sizeof(REAL4)));
      memcpy(waveform->phi->data->data ,phi.data, count*(sizeof(REAL8)));

      waveform->a->deltaT = waveform->f->deltaT = waveform->phi->deltaT
	= ppnParams->deltaT;

      waveform->a->sampleUnits    = lalStrainUnit;
      waveform->f->sampleUnits    = lalHertzUnit;
      waveform->phi->sampleUnits  = lalDimensionlessUnit;
      waveform->position = ppnParams->position;
      waveform->psi = ppnParams->psi;

      snprintf( waveform->a->name, LALNameLength,   "T1 inspiral amplitude" );
      snprintf( waveform->f->name, LALNameLength,   "T1 inspiral frequency" );
      snprintf( waveform->phi->name, LALNameLength, "T1 inspiral phase" );

      /* --- fill some output ---*/
      ppnParams->tc     = (double)(count-1) / params->tSampling ;
      ppnParams->length = count;
      ppnParams->dfdt   = ((REAL4)(waveform->f->data->data[count-1]
				   - waveform->f->data->data[count-2]))
	* ppnParams->deltaT;
      ppnParams->fStop  = params->fFinal;
      ppnParams->termCode        = GENERATEPPNINSPIRALH_EFSTOP;
      ppnParams->termDescription = GENERATEPPNINSPIRALH_MSGEFSTOP;

      ppnParams->fStart   = ppnParams->fStartIn;

  /* --- free memory --- */
  LALFree(ff.data);
  LALFree(a.data);
  LALFree(phi.data);

  DETATCHSTATUSPTR(status);
  RETURN (status);
}
/**
 * \author Creighton, T. D.
 *
 * \brief Computes a continuous waveform with frequency drift and Doppler
 * modulation from an elliptical orbital trajectory.
 *
 * This function computes a quaiperiodic waveform using the spindown and
 * orbital parameters in <tt>*params</tt>, storing the result in
 * <tt>*output</tt>.
 *
 * In the <tt>*params</tt> structure, the routine uses all the "input"
 * fields specified in \ref GenerateSpinOrbitCW_h, and sets all of the
 * "output" fields.  If <tt>params-\>f</tt>=\c NULL, no spindown
 * modulation is performed.  If <tt>params-\>oneMinusEcc</tt>\f$\notin(0,1]\f$
 * (an open orbit), or if
 * <tt>params-\>rPeriNorm</tt>\f$\times\f$<tt>params-\>angularSpeed</tt>\f$\geq1\f$
 * (faster-than-light speed at periapsis), an error is returned.
 *
 * In the <tt>*output</tt> structure, the field <tt>output-\>h</tt> is
 * ignored, but all other pointer fields must be set to \c NULL.  The
 * function will create and allocate space for <tt>output-\>a</tt>,
 * <tt>output-\>f</tt>, and <tt>output-\>phi</tt> as necessary.  The
 * <tt>output-\>shift</tt> field will remain set to \c NULL.
 *
 * ### Algorithm ###
 *
 * For elliptical orbits, we combine \eqref{eq_spinorbit-tr},
 * \eqref{eq_spinorbit-t}, and \eqref{eq_spinorbit-upsilon} to get \f$t_r\f$
 * directly as a function of the eccentric anomaly \f$E\f$:
 * \f{eqnarray}{
 * \label{eq_tr-e1}
 * t_r = t_p & + & \left(\frac{r_p \sin i}{c}\right)\sin\omega\\
 * & + & \left(\frac{P}{2\pi}\right) \left( E +
 * \left[v_p(1-e)\cos\omega - e\right]\sin E
 * + \left[v_p\sqrt{\frac{1-e}{1+e}}\sin\omega\right]
 * [\cos E - 1]\right) \;,
 * \f}
 * where \f$v_p=r_p\dot{\upsilon}_p\sin i/c\f$ is a normalized velocity at
 * periapsis and \f$P=2\pi\sqrt{(1+e)/(1-e)^3}/\dot{\upsilon}_p\f$ is the
 * period of the orbit.  For simplicity we write this as:
 * \f{equation}{
 * \label{eq_tr-e2}
 * t_r = T_p + \frac{1}{n}\left( E + A\sin E + B[\cos E - 1] \right) \;,
 * \f}
 *
 * \figure{inject_eanomaly,eps,0.23,Function to be inverted to find eccentric anomaly}
 *
 * where \f$T_p\f$ is the \e observed time of periapsis passage and
 * \f$n=2\pi/P\f$ is the mean angular speed around the orbit.  Thus the key
 * numerical procedure in this routine is to invert the expression
 * \f$x=E+A\sin E+B(\cos E - 1)\f$ to get \f$E(x)\f$.  We note that
 * \f$E(x+2n\pi)=E(x)+2n\pi\f$, so we only need to solve this expression in
 * the interval \f$[0,2\pi)\f$, sketched to the right.
 *
 * We further note that \f$A^2+B^2<1\f$, although it approaches 1 when
 * \f$e\rightarrow1\f$, or when \f$v_p\rightarrow1\f$ and either \f$e=0\f$ or
 * \f$\omega=\pi\f$.  Except in this limit, Newton-Raphson methods will
 * converge rapidly for any initial guess.  In this limit, though, the
 * slope \f$dx/dE\f$ approaches zero at the point of inflection, and an
 * initial guess or iteration landing near this point will send the next
 * iteration off to unacceptably large or small values.  However, by
 * restricting all initial guesses and iterations to the domain
 * \f$E\in[0,2\pi)\f$, one will always end up on a trajectory branch that
 * will converge uniformly.  This should converge faster than the more
 * generically robust technique of bisection. (Note: the danger with Newton's method
 * has been found to be unstable for certain binary orbital parameters. So if
 * Newton's method fails to converge, a bisection algorithm is employed.)
 *
 * In this algorithm, we start the computation with an arbitrary initial
 * guess of \f$E=0\f$, and refine it until the we get agreement to within
 * 0.01 parts in part in \f$N_\mathrm{cyc}\f$ (where \f$N_\mathrm{cyc}\f$ is the
 * larger of the number of wave cycles in an orbital period, or the
 * number of wave cycles in the entire waveform being generated), or one
 * part in \f$10^{15}\f$ (an order of magnitude off the best precision
 * possible with \c REAL8 numbers).  The latter case indicates that
 * \c REAL8 precision may fail to give accurate phasing, and one
 * should consider modeling the orbit as a set of Taylor frequency
 * coefficients \'{a} la <tt>LALGenerateTaylorCW()</tt>.  On subsequent
 * timesteps, we use the previous timestep as an initial guess, which is
 * good so long as the timesteps are much smaller than an orbital period.
 * This sequence of guesses will have to readjust itself once every orbit
 * (as \f$E\f$ jumps from \f$2\pi\f$ down to 0), but this is relatively
 * infrequent; we don't bother trying to smooth this out because the
 * additional tests would probably slow down the algorithm overall.
 *
 * Once a value of \f$E\f$ is found for a given timestep in the output
 * series, we compute the system time \f$t\f$ via \eqref{eq_spinorbit-t},
 * and use it to determine the wave phase and (non-Doppler-shifted)
 * frequency via \eqref{eq_taylorcw-freq}
 * and \eqref{eq_taylorcw-phi}.  The Doppler shift on the frequency is
 * then computed using \eqref{eq_spinorbit-upsilon}
 * and \eqref{eq_orbit-rdot}.  We use \f$\upsilon\f$ as an intermediate in
 * the Doppler shift calculations, since expressing \f$\dot{R}\f$ directly in
 * terms of \f$E\f$ results in expression of the form \f$(1-e)/(1-e\cos E)\f$,
 * which are difficult to simplify and face precision losses when
 * \f$E\sim0\f$ and \f$e\rightarrow1\f$.  By contrast, solving for \f$\upsilon\f$ is
 * numerically stable provided that the system <tt>atan2()</tt> function is
 * well-designed.
 *
 * The routine does not account for variations in special relativistic or
 * gravitational time dilation due to the elliptical orbit, nor does it
 * deal with other gravitational effects such as Shapiro delay.  To a
 * very rough approximation, the amount of phase error induced by
 * gravitational redshift goes something like \f$\Delta\phi\sim
 * fT(v/c)^2\Delta(r_p/r)\f$, where \f$f\f$ is the typical wave frequency, \f$T\f$
 * is either the length of data or the orbital period (whichever is
 * \e smaller), \f$v\f$ is the \e true (unprojected) speed at
 * periapsis, and \f$\Delta(r_p/r)\f$ is the total range swept out by the
 * quantity \f$r_p/r\f$ over the course of the observation.  Other
 * relativistic effects such as special relativistic time dilation are
 * comparable in magnitude.  We make a crude estimate of when this is
 * significant by noting that \f$v/c\gtrsim v_p\f$ but
 * \f$\Delta(r_p/r)\lesssim 2e/(1+e)\f$; we take these approximations as
 * equalities and require that \f$\Delta\phi\lesssim\pi\f$, giving:
 * \f{equation}{
 * \label{eq_relativistic-orbit}
 * f_0Tv_p^2\frac{4e}{1+e}\lesssim1 \;.
 * \f}
 * When this critereon is violated, a warning is generated.  Furthermore,
 * as noted earlier, when \f$v_p\geq1\f$ the routine will return an error, as
 * faster-than-light speeds can cause the emission and reception times to
 * be non-monotonic functions of one another.
 */
void
LALGenerateEllipticSpinOrbitCW( LALStatus             *stat,
				PulsarCoherentGW            *output,
				SpinOrbitCWParamStruc *params )
{
  UINT4 n, i;              /* number of and index over samples */
  UINT4 nSpin = 0, j;      /* number of and index over spindown terms */
  REAL8 t, dt, tPow;       /* time, interval, and t raised to a power */
  REAL8 phi0, f0, twopif0; /* initial phase, frequency, and 2*pi*f0 */
  REAL8 f, fPrev;          /* current and previous values of frequency */
  REAL4 df = 0.0;          /* maximum difference between f and fPrev */
  REAL8 phi;               /* current value of phase */
  REAL8 p, vDotAvg;        /* orbital period, and 2*pi/(period) */
  REAL8 vp;                /* projected speed at periapsis */
  REAL8 upsilon, argument; /* true anomaly, and argument of periapsis */
  REAL8 eCosOmega;         /* eccentricity * cosine of argument */
  REAL8 tPeriObs;          /* time of observed periapsis */
  REAL8 spinOff;           /* spin epoch - orbit epoch */
  REAL8 x;                 /* observed mean anomaly */
  REAL8 dx, dxMax;         /* current and target errors in x */
  REAL8 a, b;              /* constants in equation for x */
  REAL8 ecc;               /* orbital eccentricity */
  REAL8 oneMinusEcc, onePlusEcc; /* 1 - ecc and 1 + ecc */
  REAL8 e = 0.0;                 /* eccentric anomaly */
  REAL8 de = 0.0;                /* eccentric anomaly step */
  REAL8 sine = 0.0, cose = 0.0;  /* sine of e, and cosine of e minus 1 */
  REAL8 *fSpin = NULL;           /* pointer to Taylor coefficients */
  REAL4 *fData;                  /* pointer to frequency data */
  REAL8 *phiData;                /* pointer to phase data */

  INITSTATUS(stat);
  ATTATCHSTATUSPTR( stat );

  /* Make sure parameter and output structures exist. */
  ASSERT( params, stat, GENERATESPINORBITCWH_ENUL,
	  GENERATESPINORBITCWH_MSGENUL );
  ASSERT( output, stat, GENERATESPINORBITCWH_ENUL,
	  GENERATESPINORBITCWH_MSGENUL );

  /* Make sure output fields don't exist. */
  ASSERT( !( output->a ), stat, GENERATESPINORBITCWH_EOUT,
	  GENERATESPINORBITCWH_MSGEOUT );
  ASSERT( !( output->f ), stat, GENERATESPINORBITCWH_EOUT,
	  GENERATESPINORBITCWH_MSGEOUT );
  ASSERT( !( output->phi ), stat, GENERATESPINORBITCWH_EOUT,
	  GENERATESPINORBITCWH_MSGEOUT );
  ASSERT( !( output->shift ), stat, GENERATESPINORBITCWH_EOUT,
	  GENERATESPINORBITCWH_MSGEOUT );

  /* If Taylor coeficients are specified, make sure they exist. */
  if ( params->f ) {
    ASSERT( params->f->data, stat, GENERATESPINORBITCWH_ENUL,
	    GENERATESPINORBITCWH_MSGENUL );
    nSpin = params->f->length;
    fSpin = params->f->data;
  }

  /* Set up some constants (to avoid repeated calculation or
     dereferencing), and make sure they have acceptable values. */
  oneMinusEcc = params->oneMinusEcc;
  ecc = 1.0 - oneMinusEcc;
  onePlusEcc = 1.0 + ecc;
  if ( ecc < 0.0 || oneMinusEcc <= 0.0 ) {
    ABORT( stat, GENERATESPINORBITCWH_EECC,
	   GENERATESPINORBITCWH_MSGEECC );
  }
  vp = params->rPeriNorm*params->angularSpeed;
  n = params->length;
  dt = params->deltaT;
  f0 = fPrev = params->f0;
  vDotAvg = params->angularSpeed
    *sqrt( oneMinusEcc*oneMinusEcc*oneMinusEcc/onePlusEcc );
  if ( vp >= 1.0 ) {
    ABORT( stat, GENERATESPINORBITCWH_EFTL,
	   GENERATESPINORBITCWH_MSGEFTL );
  }
  if ( vp <= 0.0 || dt <= 0.0 || f0 <= 0.0 || vDotAvg <= 0.0 ||
       n == 0 ) {
    ABORT( stat, GENERATESPINORBITCWH_ESGN,
	   GENERATESPINORBITCWH_MSGESGN );
  }

  /* Set up some other constants. */
  twopif0 = f0*LAL_TWOPI;
  phi0 = params->phi0;
  argument = params->omega;
  p = LAL_TWOPI/vDotAvg;
  a = vp*oneMinusEcc*cos( argument ) + oneMinusEcc - 1.0;
  b = vp*sqrt( oneMinusEcc/( onePlusEcc ) )*sin( argument );
  eCosOmega = ecc*cos( argument );
  if ( n*dt > p )
    dxMax = 0.01/( f0*n*dt );
  else
    dxMax = 0.01/( f0*p );
  if ( dxMax < 1.0e-15 ) {
    dxMax = 1.0e-15;
#ifndef NDEBUG
    LALWarning( stat, "REAL8 arithmetic may not have sufficient"
		" precision for this orbit" );
#endif
  }
#ifndef NDEBUG
  if ( lalDebugLevel & LALWARNING ) {
    REAL8 tau = n*dt;
    if ( tau > p )
      tau = p;
    if ( f0*tau*vp*vp*ecc/onePlusEcc > 0.25 )
      LALWarning( stat, "Orbit may have significant relativistic"
		  " effects that are not included" );
  }
#endif

  /* Compute offset between time series epoch and observed periapsis,
     and betweem true periapsis and spindown reference epoch. */
  tPeriObs = (REAL8)( params->orbitEpoch.gpsSeconds -
		      params->epoch.gpsSeconds );
  tPeriObs += 1.0e-9 * (REAL8)( params->orbitEpoch.gpsNanoSeconds -
				params->epoch.gpsNanoSeconds );
  tPeriObs += params->rPeriNorm*sin( params->omega );
  spinOff = (REAL8)( params->orbitEpoch.gpsSeconds -
		     params->spinEpoch.gpsSeconds );
  spinOff += 1.0e-9 * (REAL8)( params->orbitEpoch.gpsNanoSeconds -
			       params->spinEpoch.gpsNanoSeconds );

  /* Allocate output structures. */
  if ( ( output->a = (REAL4TimeVectorSeries *)
	 LALMalloc( sizeof(REAL4TimeVectorSeries) ) ) == NULL ) {
    ABORT( stat, GENERATESPINORBITCWH_EMEM,
	   GENERATESPINORBITCWH_MSGEMEM );
  }
  memset( output->a, 0, sizeof(REAL4TimeVectorSeries) );
  if ( ( output->f = (REAL4TimeSeries *)
	 LALMalloc( sizeof(REAL4TimeSeries) ) ) == NULL ) {
    LALFree( output->a ); output->a = NULL;
    ABORT( stat, GENERATESPINORBITCWH_EMEM,
	   GENERATESPINORBITCWH_MSGEMEM );
  }
  memset( output->f, 0, sizeof(REAL4TimeSeries) );
  if ( ( output->phi = (REAL8TimeSeries *)
	 LALMalloc( sizeof(REAL8TimeSeries) ) ) == NULL ) {
    LALFree( output->a ); output->a = NULL;
    LALFree( output->f ); output->f = NULL;
    ABORT( stat, GENERATESPINORBITCWH_EMEM,
	   GENERATESPINORBITCWH_MSGEMEM );
  }
  memset( output->phi, 0, sizeof(REAL8TimeSeries) );

  /* Set output structure metadata fields. */
  output->position = params->position;
  output->psi = params->psi;
  output->a->epoch = output->f->epoch = output->phi->epoch
    = params->epoch;
  output->a->deltaT = n*params->deltaT;
  output->f->deltaT = output->phi->deltaT = params->deltaT;
  output->a->sampleUnits = lalStrainUnit;
  output->f->sampleUnits = lalHertzUnit;
  output->phi->sampleUnits = lalDimensionlessUnit;
  snprintf( output->a->name, LALNameLength, "CW amplitudes" );
  snprintf( output->f->name, LALNameLength, "CW frequency" );
  snprintf( output->phi->name, LALNameLength, "CW phase" );

  /* Allocate phase and frequency arrays. */
  LALSCreateVector( stat->statusPtr, &( output->f->data ), n );
  BEGINFAIL( stat ) {
    LALFree( output->a );   output->a = NULL;
    LALFree( output->f );   output->f = NULL;
    LALFree( output->phi ); output->phi = NULL;
  } ENDFAIL( stat );
  LALDCreateVector( stat->statusPtr, &( output->phi->data ), n );
  BEGINFAIL( stat ) {
    TRY( LALSDestroyVector( stat->statusPtr, &( output->f->data ) ),
	 stat );
    LALFree( output->a );   output->a = NULL;
    LALFree( output->f );   output->f = NULL;
    LALFree( output->phi ); output->phi = NULL;
  } ENDFAIL( stat );

  /* Allocate and fill amplitude array. */
  {
    CreateVectorSequenceIn in; /* input to create output->a */
    in.length = 2;
    in.vectorLength = 2;
    LALSCreateVectorSequence( stat->statusPtr, &(output->a->data), &in );
    BEGINFAIL( stat ) {
      TRY( LALSDestroyVector( stat->statusPtr, &( output->f->data ) ),
	   stat );
      TRY( LALDDestroyVector( stat->statusPtr, &( output->phi->data ) ),
	   stat );
      LALFree( output->a );   output->a = NULL;
      LALFree( output->f );   output->f = NULL;
      LALFree( output->phi ); output->phi = NULL;
    } ENDFAIL( stat );
    output->a->data->data[0] = output->a->data->data[2] = params->aPlus;
    output->a->data->data[1] = output->a->data->data[3] = params->aCross;
  }

  /* Fill frequency and phase arrays. */
  fData = output->f->data->data;
  phiData = output->phi->data->data;
  for ( i = 0; i < n; i++ ) {
    INT4 nOrb; /* number of orbits since the specified orbit epoch */

    /* First, find x in the range [0,2*pi]. */
    x = vDotAvg*( i*dt - tPeriObs );
    nOrb = (INT4)( x/LAL_TWOPI );
    if ( x < 0.0 )
      nOrb -= 1;
    x -= LAL_TWOPI*nOrb;

    /* Newton-Raphson iteration to find E(x). Maximum of 100 iterations. */
    INT4 maxiter = 100, iter = 0;
    while ( iter<maxiter && fabs( dx = e + a*sine + b*cose - x ) > dxMax ) {
      iter++;
      //Make a check on the step-size so we don't step too far
      de = dx/( 1.0 + a*cose + a - b*sine );
      if ( de > LAL_PI )
        de = LAL_PI;
      else if ( de < -LAL_PI )
        de = -LAL_PI;
      e -= de;

      if ( e < 0.0 )
        e = 0.0;
      else if ( e > LAL_TWOPI )
        e = LAL_TWOPI;
      sine = sin( e );
      cose = cos( e ) - 1.0;
    }
    /* Bisection algorithm from GSL if Newton's method (above) fails to converge. */
    if (iter==maxiter && fabs( dx = e + a*sine + b*cose - x ) > dxMax ) {
       //Initialize solver
       const gsl_root_fsolver_type *T = gsl_root_fsolver_bisection;
       gsl_root_fsolver *s = gsl_root_fsolver_alloc(T);
       REAL8 e_lo = 0.0, e_hi = LAL_TWOPI;
       gsl_function F;
       struct E_solver_params pars = {a, b, x};
       F.function = &gsl_E_solver;
       F.params = &pars;

       if (gsl_root_fsolver_set(s, &F, e_lo, e_hi) != 0) {
          LALFree( output->a );   output->a = NULL;
          LALFree( output->f );   output->f = NULL;
          LALFree( output->phi ); output->phi = NULL;
          ABORT( stat, -1, "GSL failed to set initial points" );
       }

       INT4 keepgoing = 1;
       INT4 success = 0;
       INT4 root_status = keepgoing;
       e = 0.0;
       iter = 0;
       while (root_status==keepgoing && iter<maxiter) {
          iter++;
          root_status = gsl_root_fsolver_iterate(s);
          if (root_status!=keepgoing && root_status!=success) {
             LALFree( output->a );   output->a = NULL;
             LALFree( output->f );   output->f = NULL;
             LALFree( output->phi ); output->phi = NULL;
             ABORT( stat, -1, "gsl_root_fsolver_iterate() failed" );
          }
          e = gsl_root_fsolver_root(s);
          sine = sin(e);
          cose = cos(e) - 1.0;
          if (fabs( dx = e + a*sine + b*cose - x ) > dxMax) root_status = keepgoing;
          else root_status = success;
       }

       if (root_status!=success) {
          LALFree( output->a );   output->a = NULL;
          LALFree( output->f );   output->f = NULL;
          LALFree( output->phi ); output->phi = NULL;
          gsl_root_fsolver_free(s);
          ABORT( stat, -1, "Could not converge using bisection algorithm" );
       }

       gsl_root_fsolver_free(s);
    }

    /* Compute source emission time, phase, and frequency. */
    phi = t = tPow =
      ( e + LAL_TWOPI*nOrb - ecc*sine )/vDotAvg + spinOff;
    f = 1.0;
    for ( j = 0; j < nSpin; j++ ) {
      f += fSpin[j]*tPow;
      phi += fSpin[j]*( tPow*=t )/( j + 2.0 );
    }

    /* Appy frequency Doppler shift. */
    upsilon = 2.0 * atan2 ( sqrt(onePlusEcc/oneMinusEcc) * sin(0.5*e), cos(0.5*e) );
    f *= f0 / ( 1.0 + vp*( cos( argument + upsilon ) + eCosOmega ) /onePlusEcc );
    phi *= twopif0;
    if ( (i > 0) && (fabs( f - fPrev ) > df) )
      df = fabs( f - fPrev );
    *(fData++) = fPrev = f;
    *(phiData++) = phi + phi0;

  } /* for i < n */

  /* Set output field and return. */
  params->dfdt = df*dt;
  DETATCHSTATUSPTR( stat );
  RETURN( stat );
}