LALSimulation  5.4.0.1-fe68b98
LALSimIMRSpinAlignedEOBHcapDerivativeOptimized.c
Go to the documentation of this file.
1 /*
2 * Copyright (C) 2010 Craig Robinson, Yi Pan
3 *
4 * This program is free software; you can redistribute it and/or modify
5 * it under the terms of the GNU General Public License as published by
6 * the Free Software Foundation; either version 2 of the License, or
7 * (at your option) any later version.
8 *
9 * This program is distributed in the hope that it will be useful,
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 * GNU General Public License for more details.
13 *
14 * You should have received a copy of the GNU General Public License
15 * along with with program; see the file COPYING. If not, write to the
16 * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
17 * MA 02110-1301 USA
18 */
19 
20 
21 /**
22  * \author Craig Robinson, Yi Pan
23  *
24  * \brief In newer versions of the EOBNR approximant, we
25  * do not have an analytic expression for the derivative of the waveform.
26  * As such, it is necessary to calculate the derivatives numerically. This
27  * function provides the means to do that for the SEOBNRv1 Hamiltonian
28  * in the spin-aligned case, i.e. equatorial orbits with four dynamical variables:
29  * r, phi, pr and pphi. It then combines energy flux and numerical derivatives of
30  * the Hamiltonian to calculate the right hand side of the ODEs given by
31  * Eqs. 10a - 10d of Pan et al. PRD 84, 124052 (2011).
32  * Since SEOBNRv1 is a spin-aligned model, the ODEs are written
33  * in the same format as ODEs of a nonspinning model.
34  *
35  */
36 
37 #ifndef LALSIMIMRSPINALIGNEDEOBHCAPDERIVATIVEOPTIMIZED_C
38 #define LALSIMIMRSPINALIGNEDEOBHCAPDERIVATIVEOPTIMIZED_C
39 
40 #include <lal/LALSimInspiral.h>
41 #include <lal/LALSimIMR.h>
42 
43 #include "LALSimIMRSpinEOB.h"
46 
47 #include <gsl/gsl_deriv.h>
48 
51 
52 /*------------------------------------------------------------------------------------------
53  *
54  * Prototypes of functions defined in this code.
55  *
56  *------------------------------------------------------------------------------------------
57  */
58 
60  double t,
61  const REAL8 values[],
62  REAL8 dvalues[],
63  void *funcParams
64  );
65 /*------------------------------------------------------------------------------------------
66  *
67  * Defintions of functions.
68  *
69  *------------------------------------------------------------------------------------------
70  */
71 
72 /**
73  * Function to calculate R.H.S. of the ODEs, given dyanmical variables,
74  * their derivatives and EOB parameters. Since SEOBNRv1 spin Hamiltonian
75  * was implemented for Cartesean coordinates while dynamical evolution was
76  * implemented in polar coordinates, we need to perform a transform.
77  * This is done in a particular transform in which
78  * x = r, y = z = 0, px = pr, py = pphi/r, pz = 0, and
79  * omega = v/r = (dy/dt)/r = (dH/dpy)/r, dr/dt = dx/dt = dH/dpx, etc.
80  */
82  double UNUSED t, /**< UNUSED */
83  const REAL8 values[], /**< dynamical varables */
84  REAL8 dvalues[], /**< time derivative of dynamical variables */
85  void *funcParams /**< EOB parameters */
86  )
87 {
88  static const INT4 lMax = 8;
89 
91 
92  /* Since we take numerical derivatives wrt dynamical variables */
93  /* but we want them wrt time, we use this temporary vector in */
94  /* the conversion */
95  REAL8 tmpDValues[6];
96 
97  /* Cartesian values for calculating the Hamiltonian */
98  REAL8 cartValues[6];
99 
100  REAL8 H; //Hamiltonian
101  REAL8 flux;
102 
103  UINT4 SpinAlignedEOBversion;
104 
105  REAL8Vector rVec, pVec;
106  REAL8 rData[3], pData[3];
107 
108  /* We need r, phi, pr, pPhi to calculate the flux */
109  REAL8 r;
110  REAL8Vector polarDynamics;
111  REAL8 polData[4];
112 
113  REAL8 mass1, mass2, eta;
114 
115  /* Spins */
116  REAL8Vector *s1Vec = NULL;
117  REAL8Vector *s2Vec = NULL;
118  REAL8Vector *sKerr = NULL;
119  REAL8Vector *sStar = NULL;
120 
121  REAL8 a;
122 
123  REAL8 omega;
124 
125  /* EOB potential functions */
126  REAL8 DeltaT, DeltaR;
127  REAL8 csi;
128 
129  /* The error in a derivative as measured by GSL */
130 
131  /* Declare NQC coefficients */
132  EOBNonQCCoeffs *nqcCoeffs = NULL;
133 
134  /* Set up pointers for GSL */
135  params.values = cartValues;
136  params.params = (SpinEOBParams *)funcParams;
137  nqcCoeffs = params.params->nqcCoeffs;
138 
139  s1Vec = params.params->s1Vec;
140  s2Vec = params.params->s2Vec;
141  sKerr = params.params->sigmaKerr;
142  sStar = params.params->sigmaStar;
143 
144  mass1 = params.params->eobParams->m1;
145  mass2 = params.params->eobParams->m2;
146  eta = params.params->eobParams->eta;
147 
148  SpinAlignedEOBversion = params.params->seobCoeffs->SpinAlignedEOBversion;
149 
150  r = values[0];
151 
152  /* Since this is spin aligned, I make the assumption */
153  /* that the spin vector is along the z-axis. */
154  a = sKerr->data[2];
155 
156  /* Calculate the potential functions and the tortoise coordinate factor csi,
157  given by Eq. 28 of Pan et al. PRD 81, 084041 (2010) */
158  DeltaT = XLALSimIMRSpinEOBHamiltonianDeltaT( params.params->seobCoeffs, r, eta, a );
159  DeltaR = XLALSimIMRSpinEOBHamiltonianDeltaR( params.params->seobCoeffs, r, eta, a );
160  csi = sqrt( DeltaT * DeltaR ) / (r*r + a*a);
161  //printf("DeltaT = %.16e, DeltaR = %.16e, a = %.16e\n",DeltaT,DeltaR,a);
162  //printf( "csi in derivatives function = %.16e\n", csi );
163 
164  /* Populate the Cartesian values vector, using polar coordinate values */
165  /* We can assume phi is zero wlog */
166  memset( cartValues, 0, sizeof( cartValues ) );
167  cartValues[0] = values[0];
168  cartValues[3] = values[2];
169  cartValues[4] = values[3] / values[0];
170 
171  /* Now calculate derivatives w.r.t. each Cartesian variable */
172  /* OPTIMIZATION NOTE:
173  * The original function computed derivatives of the Hamiltonian with respect to 6 variables (three coordinate & three momenta).
174  * In SEOBNRv2, only three are nonzero! */
175  GSLSpinAlignedHamiltonianWrapper_derivs_allatonce( tmpDValues, cartValues, &params );
176 
177  /* Calculate the Cartesian vectors rVec and pVec */
178  polarDynamics.length = 4;
179  polarDynamics.data = polData;
180 
181  memcpy( polData, values, sizeof( polData ) );
182 
183  rVec.length = pVec.length = 3;
184  rVec.data = rData;
185  pVec.data = pData;
186 
187  memset( rData, 0, sizeof(rData) );
188  memset( pData, 0, sizeof(pData) );
189 
190  rData[0] = values[0];
191  pData[0] = values[2];
192  pData[1] = values[3] / values[0];
193  /* Calculate Hamiltonian using Cartesian vectors rVec and pVec */
194  H = XLALSimIMRSpinEOBHamiltonianOptimized( eta, &rVec, &pVec, s1Vec, s2Vec, sKerr, sStar, params.params->tortoise, params.params->seobCoeffs );
195 
196  //printf( "csi = %.16e, ham = %.16e ( tortoise = %d)\n", csi, H, params.params->tortoise );
197  //exit(1);
198  //if ( values[0] > 1.3 && values[0] < 3.9 ) printf( "r = %e\n", values[0] );
199  //if ( values[0] > 1.3 && values[0] < 3.9 ) printf( "Hamiltonian = %e\n", H );
200  /* OPTIMIZATION NOTE:
201  * Why are we multiplying and then dividing by the total mass two lines later? */
202  H = H * (mass1 + mass2);
203 
204 
205  /*if ( values[0] > 1.3 && values[0] < 3.9 ) printf( "Cartesian derivatives:\n%f %f %f %f %f %f\n",
206  tmpDValues[3], tmpDValues[4], tmpDValues[5], -tmpDValues[0], -tmpDValues[1], -tmpDValues[2] );*/
207 
208  /* Now calculate omega, and hence the flux */
209  omega = tmpDValues[4] / r;
210  flux = XLALInspiralSpinFactorizedFluxOptimized( &polarDynamics, nqcCoeffs, omega, params.params, H/(mass1+mass2), lMax, SpinAlignedEOBversion );
211 
212  /* Looking at the non-spinning model, I think we need to divide the flux by eta */
213  flux = flux / eta;
214 
215  //printf( "Flux in derivatives function = %.16e\n", flux );
216 
217  /* Now we can calculate the final (spherical) derivatives */
218  /* csi is needed because we use the tortoise co-ordinate */
219  /* Right hand side of Eqs. 10a - 10d of Pan et al. PRD 84, 124052 (2011) */
220  dvalues[0] = csi * tmpDValues[3];
221  dvalues[1] = omega;
222  /* Note: in this special coordinate setting, namely y = z = 0, dpr/dt = dpx/dt + dy/dt * py/r, where py = pphi/r */
223  dvalues[2] = - tmpDValues[0] + tmpDValues[4] * values[3] / (r*r);
224  dvalues[2] = dvalues[2] * csi - ( values[2] / values[3] ) * flux / omega;
225  dvalues[3] = - flux / omega;
226 
227  //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] );
228 
229  //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] );
230 
231  if ( isnan( dvalues[0] ) || isnan( dvalues[1] ) || isnan( dvalues[2] ) || isnan( dvalues[3] ) )
232  {
233  //printf( "Deriv is nan: %e %e %e %e\n", dvalues[0], dvalues[1], dvalues[2], dvalues[3] );
234  return 1;
235  }
236 
237  return XLAL_SUCCESS;
238 }
239 
240 #endif /* LALSIMIMRSPINALIGNEDEOBHCAPDERIVATIVEOPTIMIZED_C */
static int XLALSpinAlignedHcapDerivativeOptimized(double t, const REAL8 values[], REAL8 dvalues[], void *funcParams)
static REAL8 XLALInspiralSpinFactorizedFluxOptimized(REAL8Vector *values, EOBNonQCCoeffs *nqcCoeffs, const REAL8 omega, SpinEOBParams *ak, const REAL8 H, const INT4 lMax, const UINT4 SpinAlignedEOBversion)
This function calculates the spin factorized-resummed GW energy flux for given dynamical variables.
static REAL8 XLALSimIMRSpinEOBHamiltonianDeltaR(SpinEOBHCoeffs *coeffs, const REAL8 r, const REAL8 eta, const REAL8 a)
This function calculates the DeltaR potential function in the spin EOB Hamiltonian.
static REAL8 XLALSimIMRSpinEOBHamiltonianDeltaT(SpinEOBHCoeffs *coeffs, const REAL8 r, const REAL8 eta, const REAL8 a)
This function calculates the function which appears in the spinning EOB potential function.
UNUSED REAL8 XLALSimIMRSpinEOBHamiltonianOptimized(const REAL8 eta, REAL8Vector *restrict x, REAL8Vector *restrict p, REAL8Vector *restrict s1Vec, REAL8Vector *restrict s2Vec, REAL8Vector *restrict sigmaKerr, REAL8Vector *restrict sigmaStar, INT4 tortoise, SpinEOBHCoeffs *coeffs)
Function to calculate the value of the spinning Hamiltonian for given values of the dynamical variabl...
static REAL8 GSLSpinAlignedHamiltonianWrapper_derivs_allatonce(REAL8 output[6], const REAL8 input[6], void *params)
const double H
const double csi
double REAL8
uint32_t UINT4
int32_t INT4
static const INT4 r
static const INT4 a
XLAL_SUCCESS
The coefficients which are used in calculating the non-quasicircular correction to the EOBNRv2 model.
REAL8 * data
Parameters for the spinning EOB model.
Definition: burst.c:245
LALDict * params
Definition: inspiral.c:248