int MultidimensionalRootFinder(const int dimension, gsl_multiroot_function *f, gsl_vector *initial_guess, double abs_error, double rel_error, int max_iterations, gsl_vector *results) { int status; size_t iter = 0; const gsl_multiroot_fsolver_type * solver_type = gsl_multiroot_fsolver_broyden; gsl_multiroot_fsolver * solver = gsl_multiroot_fsolver_alloc(solver_type, dimension); gsl_multiroot_fsolver_set(solver, f, initial_guess); do { iter++; status = gsl_multiroot_fsolver_iterate(solver); if (status == GSL_EBADFUNC){ printf("TwodimensionalRootFinder: Error: Infinity or division by zero.\n"); abort(); } else if (status == GSL_ENOPROG){ printf("TwodimensionalRootFinder: Error: Solver is stuck. Try a different initial guess.\n"); abort(); } // Check if the root is good enough: // tests for the convergence of the sequence by comparing the last step dx with the // absolute error epsabs and relative error epsrel to the current position x. The test // returns GSL_SUCCESS if the following condition is achieved, // // |dx_i| < epsabs + epsrel |x_i| gsl_vector * x = gsl_multiroot_fsolver_root(solver); // current root gsl_vector * dx = gsl_multiroot_fsolver_dx(solver); // last step status = gsl_multiroot_test_delta(dx, x, abs_error, rel_error); } while (status == GSL_CONTINUE && iter < max_iterations); // Save results in return variables gsl_vector_memcpy(results, gsl_multiroot_fsolver_root(solver)); // Free vectors gsl_multiroot_fsolver_free(solver); 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; }