precession.eccentricity: package documentation
- eccentricity.affine(vec, low, up)[source]
Perform an affine trasformation to scale a vector given upper and lower limits: rescaled = ( vec - low ) / (up - low). If these are outside the domain of the vector, the result will be in [0,1].
- Parameters:
- vec: array
Vector in Cartesian coomponents.
- low: float
Lower limit.
- up: float
Upper limit.
- Returns:
- rescaled: float, optional
Rescaled vector.
Examples
rescaled = precession.affine(vec,low,up)
- eccentricity.angleresonances_endpoint(q, chi1, chi2, chieff)[source]
Small-separation limit of the spin-orbit resonances.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- chieff: float
Effective spin.
- Returns:
- deltaphiatmax: float
Value of the angle deltaphi at the resonance that maximizes kappa.
- deltaphiatmin: float
Value of the angle deltaphi at the resonance that minimizes kappa.
- theta1atmax: float
Value of the angle theta1 at the resonance that maximizes kappa.
- theta1atmin: float
Value of the angle theta1 at the resonance that minimizes kappa.
- theta2atmax: float
Value of the angle theta2 at the resonance that maximizes kappa.
- theta2atmin: float
Value of the angle theta2 at the resonance that minimizes kappa.
Examples
theta1atmin,theta2atmin,deltaphiatmin,theta1atmax,theta2atmax,deltaphiatmax = precession.angleresonances_endpoint(q,chi1,chi2,chieff)
- eccentricity.angles_to_Jframe(theta1, theta2, deltaphi, a, e, q, chi1, chi2)[source]
Convert the angles (theta1,theta2,deltaphi) to angular momentum vectors (L,S1,S2) in the frame aligned with the total angular momentum. In particular, we set Jx=Jy=Ly=0.
- Parameters:
- theta1: float
Angle between orbital angular momentum and primary spin.
- theta2: float
Angle between orbital angular momentum and secondary spin.
- deltaphi: float
Angle between the projections of the two spins onto the orbital plane.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- Lvec: array
Cartesian vector of the orbital angular momentum.
- S1vec: array
Cartesian vector of the primary spin.
- S2vec: array
Cartesian vector of the secondary spin.
Examples
Lvec,S1vec,S2vec = precession.angles_to_Jframe(theta1,theta2,deltaphi,a,e,q,chi1,chi2)
- eccentricity.angles_to_Lframe(theta1, theta2, deltaphi, a, e, q, chi1, chi2)[source]
Convert the angles (theta1,theta2,deltaphi) to angular momentum vectors (L,S1,S2) in the frame aligned with the orbital angular momentum. In particular, we set Lx=Ly=S1y=0.
- Parameters:
- theta1: float
Angle between orbital angular momentum and primary spin.
- theta2: float
Angle between orbital angular momentum and secondary spin.
- deltaphi: float
Angle between the projections of the two spins onto the orbital plane.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- Lvec: array
Cartesian vector of the orbital angular momentum.
- S1vec: array
Cartesian vector of the primary spin.
- S2vec: array
Cartesian vector of the secondary spin.
Examples
Lvec,S1vec,S2vec = precession.angles_to_Lframe(theta1,theta2,deltaphi,a,e,q,chi1,chi2)
- eccentricity.angles_to_conserved(theta1, theta2, deltaphi, a, e, q, chi1, chi2, full_output=False)[source]
Convert angles (theta1,theta2,deltaphi) into conserved quantities (deltachi,kappa,chieff).
- Parameters:
- theta1: float
Angle between orbital angular momentum and primary spin.
- theta2: float
Angle between orbital angular momentum and secondary spin.
- deltaphi: float
Angle between the projections of the two spins onto the orbital plane.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- full_output: boolean, optional (default: False)
Return additional outputs.
- Returns:
- chieff: float
Effective spin.
- cyclesign: integer, optional
Sign (either +1 or -1) to cover the two halves of a precesion cycle.
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
Examples
deltachi,kappa,chieff = precession.angles_to_conserved(theta1,theta2,deltaphi,a,e,q,chi1,chi2)
deltachi,kappa,chieff,cyclesign = precession.angles_to_conserved(theta1,theta2,deltaphi,a,e,q,chi1,chi2,full_output=True)
- eccentricity.anglesresonances(a, e, chieff, q, chi1, chi2)[source]
Compute the values of the angles corresponding to the two spin-orbit resonances.
- Parameters:
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- deltaphiatmax: float
Value of the angle deltaphi at the resonance that maximizes kappa.
- deltaphiatmin: float
Value of the angle deltaphi at the resonance that minimizes kappa.
- theta1atmax: float
Value of the angle theta1 at the resonance that maximizes kappa.
- theta1atmin: float
Value of the angle theta1 at the resonance that minimizes kappa.
- theta2atmax: float
Value of the angle theta2 at the resonance that maximizes kappa.
- theta2atmin: float
Value of the angle theta2 at the resonance that minimizes kappa.
Examples
theta1atmin,theta2atmin,deltaphiatmin,theta1atmax,theta2atmax,deltaphiatmax = precession.anglesresonances(r,chieff,q,chi1,chi2)
- eccentricity.chiefflimits(q, chi1, chi2)[source]
Limits on the effective spin based only on its definition.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- chieffmax: float
Maximum value of the effective spin.
- chieffmin: float
Minimum value of the effective spin.
Examples
chieffmin,chieffmax = precession.chiefflimits(q,chi1,chi2)
- eccentricity.chip_terms(theta1, theta2, q, chi1, chi2)[source]
Compute the two terms entering the effective precessing spin chip.
- Parameters:
- theta1: float
Angle between orbital angular momentum and primary spin.
- theta2: float
Angle between orbital angular momentum and secondary spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- chipterm1: float
Term in effective precessing spin.
- chipterm2: float
Term in effective precessing spin.
Examples
chipterm1,chipterm2 = precession.chip_terms(theta1,theta2,q,chi1,chi2)
- eccentricity.conserved_to_Jframe(deltachi, kappa, a, e, chieff, q, chi1, chi2, cyclesign=1)[source]
Convert the conserved quanties (deltachi,kappa,chieff) to angular momentum vectors (L,S1,S2) in the frame aligned with the total angular momentum. In particular, we set Jx=Jy=Ly=0.
- Parameters:
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- cyclesign: integer, optional (default: +1)
Sign (either +1 or -1) to cover the two halves of a precesion cycle.
- Returns:
- Lvec: array
Cartesian vector of the orbital angular momentum.
- S1vec: array
Cartesian vector of the primary spin.
- S2vec: array
Cartesian vector of the secondary spin.
Examples
Lvec,S1vec,S2vec = precession.conserved_to_Jframe(deltachi,kappa,a,e,chieff,q,chi1,chi2,cyclesign=+1)
- eccentricity.conserved_to_Lframe(deltachi, kappa, a, e, chieff, q, chi1, chi2, cyclesign=1)[source]
Convert the conserved quanties (deltachi,kappa,chieff) to angular momentum vectors (L,S1,S2) in the frame aligned with the orbital angular momentum. In particular, we set Lx=Ly=S1y=0.
- Parameters:
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- cyclesign: integer, optional (default: +1)
Sign (either +1 or -1) to cover the two halves of a precesion cycle.
- Returns:
- Lvec: array
Cartesian vector of the orbital angular momentum.
- S1vec: array
Cartesian vector of the primary spin.
- S2vec: array
Cartesian vector of the secondary spin.
Examples
Lvec,S1vec,S2vec = precession.conserved_to_Lframe(deltachi,kappa,a,e,chieff,q,chi1,chi2,cyclesign=+1)
- eccentricity.conserved_to_angles(deltachi, kappa, a, e, chieff, q, chi1, chi2, cyclesign=1)[source]
Convert conserved quantities (deltachi,kappa,chieff) into angles (theta1,theta2,deltaphi).
- Parameters:
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- cyclesign: integer, optional (default: +1)
Sign (either +1 or -1) to cover the two halves of a precesion cycle.
- Returns:
- deltaphi: float
Angle between the projections of the two spins onto the orbital plane.
- theta1: float
Angle between orbital angular momentum and primary spin.
- theta2: float
Angle between orbital angular momentum and secondary spin.
Examples
theta1,theta2,deltaphi = precession.conserved_to_angles(deltachi,kappa,a,e,chieff,q,chi1,chi2,cyclesign=+1)
- eccentricity.dchidt2_RHS(deltachi, kappa, a, e, chieff, q, chi1, chi2, precomputedroots=None, donotnormalize=False)[source]
Right-hand side of the (ddeltachi/dt)^2 equation. Setting donotnormalize=True is equivalent to mathcalA=1.
- Parameters:
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- donotnormalize: boolean, optional (default: False)
If True omit the numerical prefactor.
- Returns:
- dchidt2: float
Squared time derivative of the weighted spin difference.
Examples
dchidt2 = precession.dchidt2_RHS(r,chieff,q)
- eccentricity.ddchidt_prefactor(a, e, chieff, q)[source]
etamerical prefactor to the ddeltachi/dt derivative.
- Parameters:
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1.
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- Returns:
- mathcalA: float
Prefactor in the ddeltachi/dt equation.
Examples
mathcalA = precession.ddchidt_prefactor(a,e,chieff,q)
- eccentricity.deltachicubic_coefficients(kappa, u, chieff, q, chi1, chi2)[source]
Coefficients of the cubic equation in deltachi for its time evolution.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- u: float
Compactified separation 1/(2L).
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- coeff0: float
Coefficient to the x^0 term in polynomial.
- coeff1: float
Coefficient to the x^1 term in polynomial.
- coeff2: float
Coefficient to the x^2 term in polynomial.
- coeff3: float
Coefficient to the x^3 term in polynomial.
Examples
coeff3,coeff2,coeff1,coeff0 = precession.deltachicubic_coefficients(kappa,u,chieff,q,chi1,chi2)
- eccentricity.deltachicubic_rescaled_coefficients(kappa, u, chieff, q, chi1, chi2, precomputedcoefficients=None)[source]
Rescaled coefficients of the cubic equation in deltachi for its time evolution. This is necessary to avoid dividing by (1-q).
- Parameters:
- kappa: float
Asymptotic angular momentum.
- u: float
Compactified separation 1/(2L).
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedcoefficients: array, optional (default: None)
Pre-computed output of deltachicubic_coefficients for computational efficiency.
- Returns:
- coeff0: float
Coefficient to the x^0 term in polynomial.
- coeff1: float
Coefficient to the x^1 term in polynomial.
- coeff2: float
Coefficient to the x^2 term in polynomial.
- coeff3: float
Coefficient to the x^3 term in polynomial.
Examples
coeff3,coeff2,coeff1,coeff0 = precession.deltachicubic_rescaled_coefficients(kappa,u,chieff,q,chi1,chi2)
coeff3,coeff2,coeff1,coeff0 = precession.deltachicubic_rescaled_coefficients(kappa,u,chieff,q,chi1,chi2,precomputedcoefficients=coeffs)
- eccentricity.deltachilimits(kappa=None, a=None, e=None, chieff=None, q=None, chi1=None, chi2=None, precomputedroots=None)[source]
- Limits on the weighted spin difference. The contraints considered depend on the inputs provided.
If q, chi1, chi2 are provided, returns the geometrical limits.
If chieff, q, chi1, and chi2 are provided, returns the joint ‘rectagle’ limits.
If kappa, r, chieff, q, chi1, and chi2 are provided, returns the solution of the cubic that sets the precession-timescale evolution.
- Parameters:
- kappa: float, optional (default: None)
Asymptotic angular momentum.
- a: float, optional (default: None)
Semi-major axis.
- e: float, optional (default: None)
Eccentricity: 0<=e<1
- chieff: float, optional (default: None)
Effective spin.
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- chi1: float, optional (default: None)
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float, optional (default: None)
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- deltachimax: float
Maximum value of the weighted spin difference.
- deltachimin: float
Minimum value of the weighted spin difference.
Examples
deltachimin,deltachimax = precession.deltachilimits(q=q,chi1=chi1,chi2=chi2)
deltachimin,deltachimax = precession.deltachilimits(chieff=chieff,q=q,chi1=chi1,chi2=chi2)
deltachimin,deltachimax = precession.deltachilimits(kappa=kappa,a=a,e=e,chieff=chieff,q=q,chi1=chi1,chi2=chi2)
- eccentricity.deltachilimits_definition(q, chi1, chi2)[source]
Limits on the weighted spin difference based only on its definition.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- deltachimax: float
Maximum value of the effective spin chieff.
- deltachimin: float
Minimum value of the effective spin chieff.
Examples
deltachimin,deltachimax = precession.deltachilimits_definition(q,chi1,chi2)
- eccentricity.deltachilimits_plusminus(kappa, a, e, chieff, q, chi1, chi2, precomputedroots=None)[source]
Limits on the weighted spin difference for given (kappa, r, chieff, q, chi1, chi2). These are two of the solutions of the underlying cubic polynomial.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- deltachimax: float
Maximum value of the weighted spin difference.
- deltachimin: float
Minimum value of the weighted spin difference.
Examples
deltachimin,deltachimax = precession.deltachilimits_plusminus(kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.deltachilimits_rectangle(chieff, q, chi1, chi2)[source]
Limits on the weighted spin difference for given (chieff, q, chi1, chi2). This is a rectangle in the (chieff,deltachi) plane.
- Parameters:
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- deltachimax: float
Maximum value of the weighted spin difference.
- deltachimin: float
Minimum value of the weighted spin difference.
Examples
deltachimin,deltachimax = precession.deltachilimits_rectangle(chieff,q,chi1,chi2)
- eccentricity.deltachioft(t, kappa, a, e, chieff, q, chi1, chi2, precomputedroots=None)[source]
Evolution of deltachi on the precessional timescale (without radiation reaction). The broadcasting rules for this function are more general than those of the rest of the code. The variable t is allowed to have shapes (N,M) while all the other variables have shape (N,). This is useful to sample M precession configuration for each of the N binaries specified as inputs.
- Parameters:
- t: float
Time.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- deltachi: float
Weighted spin difference.
Examples
deltachi = precession.deltachioft(t,kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.deltachirescaling(deltachitilde, kappa, a, e, chieff, q, chi1, chi2, precomputedroots=None)[source]
Compute deltachi from the rescaled parameter 0<=deltachitilde<=1.
- Parameters:
- deltachitilde: float
Rescaled version of the weighted spin difference.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- deltachi: float
Weighted spin difference.
Examples
deltachi = precession.deltachirescaling(deltachitilde,kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.deltachiresonance(kappa=None, a=None, e=None, u=None, chieff=None, q=None, chi1=None, chi2=None)[source]
Assuming that the inputs correspond to a spin-orbit resonance, find the corresponding value of deltachi. There will be two roots that are coincident if not for numerical errors: for concreteness, return the mean of the real part. This function is dangerous: we do not check that the input is a resonance, that’s up to the user (so use at your own risk). Valid inputs are either (kappa,a,e,chieff,q,chi1,chi2) or (kappa,u,chieff,q,chi1,chi2).
- Parameters:
- kappa: float, optional (default: None)
Asymptotic angular momentum.
- a: float, optional (default: None)
Semi-major axis.
- e: float, optional (default: None)
Eccentricity: 0<=e<1
- u: float, optional (default: None)
Compactified separation 1/(2L).
- chieff: float, optional (default: None)
Effective spin.
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- chi1: float, optional (default: None)
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float, optional (default: None)
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- deltachi: float
Weighted spin difference.
Examples
deltachi = precession.deltachiresonance(kappa=kappa,a=a,e=e,chieff=chieff,q=q,chi1=chi1,chi2=chi2)
deltachi = precession.deltachiresonance(kappa=kappa,u=u,chieff=chieff,q=q,chi1=chi1,chi2=chi2)
- eccentricity.deltachiroots(kappa, u, chieff, q, chi1, chi2, full_output=True, precomputedroots=None)[source]
Roots of the cubic equation in deltachi that describes the dynamics on the precession timescale.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- u: float
Compactified separation 1/(2L).
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- full_output: boolean, optional (default: True)
Return additional outputs.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- deltachi3: float, optional
Spurious root of the deltachi evolution.
- deltachiminus: float
Lowest physical root of the deltachi evolution.
- deltachiplus: float
Lowest physical root of the deltachi evolution.
Examples
deltachiminus,deltachiplus,deltachi3 = precession.deltachiroots(kappa,u,chieff,q,chi1,chi2)
deltachiminus,deltachiplus,deltachi3 = precession.deltachiroots(kappa,u,chieff,q,chi1,chi2,precomputedroots=roots)
deltachiminus,deltachiplus = precession.deltachiroots(kappa,u,chieff,q,chi1,chi2,full_output=False)
- eccentricity.deltachisampling(kappa, a, e, chieff, q, chi1, chi2, N=1, precomputedroots=None)[source]
Sample N values of deltachi at fixed separation accoring to its PN-weighted distribution function. Can only be used to sample the same number of configuration for each binary. If the inputs have shape (M,) the output will have shape
(M,N) if M>1 and N>1;
(M,) if N=1;
(N,) if M=1.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- N: integer, optional (default: 1)
Number of samples.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- deltachi: float
Weighted spin difference.
Examples
deltachi = precession.deltachisampling(kappa,a,e,chieff,q,chi1,chi2,N=1)
- eccentricity.deltachitildeav(m, tol=1e-07)[source]
Precession average of the rescaled weighted spin difference.
- Parameters:
- m: float
Parameter of elliptic function(s).
- tol: float, optional (default: 1e-7)
Numerical tolerance, see source code for details.
- Returns:
- coeff: float
Coefficient.
Examples
coeff = precession.deltachitildeav(m,tol=1e-7)
- eccentricity.deltachitildeav2(m, tol=1e-07)[source]
Precession average of the rescaled weighted spin difference squared.
- Parameters:
- m: float
Parameter of elliptic function(s).
- tol: float, optional (default: 1e-7)
Numerical tolerance, see source code for details.
- Returns:
- coeff: float
Coefficient.
Examples
coeff = precession.deltachitildeav2(m,tol=1e-7)
- eccentricity.dot_nested(x, y)[source]
Dot product between 2D arrays along last axis.
- Parameters:
- xarray
Input array.
- yarray
Input array.
- Returns:
- zarray
Dot product array.
Examples
z = dot_nested(x, y)
- eccentricity.ellippi(n, phi, m)[source]
Incomplete elliptic integral of the third kind. This is reconstructed using scipy’s implementation of Carlson’s R integrals (arxiv:math/9409227).
- Parameters:
- n: foat
Characheristic of the elliptic integral.
- phi: float
Amplitude of the elliptic integral.
- m: float
Parameter of the elliptic integral
- Returns:
- piintegral: float
Incomplete elliptic integral of the third kind
Examples
piintegral = precession.ellippi(n, phi, m)
- eccentricity.elliptic_parameter(kappa, u, chieff, q, chi1, chi2, precomputedroots=None)[source]
Parameter m entering elliptic functions for the evolution of deltachi.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- u: float
Compactified separation 1/(2L).
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- m: float
Parameter of elliptic function(s).
Examples
m = precession.elliptic_parameter(kappa,u,chieff,q,chi1,chi2)
- eccentricity.eval_J(theta1=None, theta2=None, deltaphi=None, kappa=None, a=None, e=None, q=None, chi1=None, chi2=None)[source]
Magnitude of the total angular momentum. Provide either (theta1,theta2,deltaphi,a,e,q,chi1,chhi2) or (kappa,a,e,q).
- Parameters:
- theta1: float, optional (default: None)
Angle between orbital angular momentum and primary spin.
- theta2: float, optional (default: None)
Angle between orbital angular momentum and secondary spin.
- deltaphi: float, optional (default: None)
Angle between the projections of the two spins onto the orbital plane.
- kappa: float, optional (default: None)
Asymptotic angular momentum.
- a: float, optional (default: None)
Semi-major axis.
- e: float, optional (default: None)
Eccentricity: 0<=e<1
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- chi1: float, optional (default: None)
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float, optional (default: None)
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- J: float
Magnitude of the total angular momentum.
Examples
J = precession.eval_J(theta1=theta1,theta2=theta2,deltaphi=deltaphi,a=a,e=e,q=q,chi1=chi1,chi2=chi2)
J = precession.eval_J(kappa=kappa,a=a,e=e,q=q,chi1=chi1,chi2=chi2)
- eccentricity.eval_L(a, e, q)[source]
Newtonian angular momentum of the binary.
- Parameters:
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- q: float
Mass ratio: 0<=q<=1.
- Returns:
- L: float
Magnitude of the Newtonian orbital angular momentum.
Examples
L = precession.eval_L(r,q)
- eccentricity.eval_OmegaL(deltachi, kappa, a, e, chieff, q, chi1, chi2)[source]
Evaluate the precession frequency OmegaL along the precession cycle.
- Parameters:
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- OmegaL: float
Precession frequency of L about J.
Examples
OmegaL = precession.eval_OmegaL(deltachi,kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.eval_S(theta1=None, theta2=None, deltaphi=None, deltachi=None, kappa=None, a=None, e=None, chieff=None, q=None, chi1=None, chi2=None)[source]
Magnitude of the total spin. Valid inputs are either (theta1, theta2, deltaphi, q, chi1, chi2) or (deltachi, kappa, r, chieff, q).
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- Returns:
- S: float
Magnitude of the total spin.
Examples
S = precession.eval_S(theta1=theta1,theta2=theta2,deltaphi=deltaphi,q=q,chi1=chi1,chi2=chi2)
S = precession.eval_S(deltachi=deltachi,kappa=kappa,a=a,e=e,chieff=chieff,q=q)
- eccentricity.eval_S1(q, chi1)[source]
Spin angular momentum of the heavier black hole.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- Returns:
- S1: float
Magnitude of the primary spin.
Examples
S1 = precession.eval_S1(q,chi1)
- eccentricity.eval_S2(q, chi2)[source]
Spin angular momentum of the lighter black hole.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- S2: float
Magnitude of the secondary spin.
Examples
S2 = precession.eval_S2(q,chi2)
- eccentricity.eval_a(L=None, u=None, uc=None, e=None, q=None)[source]
Semi-major axis of the binary. Valid inputs are either (L,e,q) or (u,e,q) or (uc,q).
- Parameters:
- L: float, optional (default: None)
Magnitude of the Newtonian orbital angular momentum.
- u: float, optional (default: None)
Compactified separation 1/(2L).
- uc: float, optional (default: None)
Circular compactified separation uc=u(e=0).
- e: float (default: 0)
Binary eccentricity: 0<=e<1.
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- Returns:
- a: float
Binary semi-major axis.
- eccentricity.eval_alpha(kappa, a, e, chieff, q, chi1, chi2, precomputedroots=None)[source]
Evaluate the precession angle spanned by L during an entire cycle.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- alpha: float
Azimuthal angle spanned by L about J during an entire cycle.
Examples
alpha = precession.eval_alpha(kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.eval_bracket_omega(kappa, a, e, chieff, q, chi1, chi2, precomputedroots=None)[source]
Precession average of the precession frequency deltachi as it oscillates from deltachi- to deltachi+ back to deltachi-. This of the five phenomenological parameters from arxiv:2103.03894.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- bracket_omega: float
Precession-averaged precession frequency.
Examples
bracket_omega = precession.eval_bracket_omega(kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.eval_bracket_theta(kappa, a, e, chieff, q, chi1, chi2, **kwargs)[source]
Precession average of precession amplitude of deltachi as it oscillates from deltachi- to deltachi+ back to deltachi-. This of the five phenomenological parameters from arxiv:2103.03894. Additional keywords arguments are passed to precession_average.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- **kwargs: unpacked dictionary, optional
Additional keyword arguments.
- Returns:
- bracket_theta: float
Precession-averaged precession amplitude.
Examples
bracket_theta = precession.eval_bracket_theta(kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.eval_chi1(q, S1)[source]
Dimensionless spin of the heavier black hole.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- S1: float
Magnitude of the primary spin.
- Returns:
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
Examples
chi1 = precession.eval_chi1(q,S1)
- eccentricity.eval_chi2(q, S2)[source]
Dimensionless spin of the lighter black hole.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- S2: float
Magnitude of the secondary spin.
- Returns:
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
Examples
chi2 = precession.eval_chi2(q,S2)
- eccentricity.eval_chieff(theta1, theta2, q, chi1, chi2)[source]
Eftective spin.
- Parameters:
- theta1: float
Angle between orbital angular momentum and primary spin.
- theta2: float
Angle between orbital angular momentum and secondary spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- chieff: float
Effective spin.
Examples
chieff = precession.eval_chieff(theta1,theta2,q,chi1,chi2)
- eccentricity.eval_chip(theta1=None, theta2=None, deltaphi=None, deltachi=None, kappa=None, a=None, e=None, chieff=None, q=None, chi1=None, chi2=None, which='averaged', **kwargs)[source]
- Compute the effective precessing spin. The keyword which one of the following definitions:
heuristic, as in Schmidt et al 2015. Required inputs are either
generalized, retail all precession-timescale variations,
averaged (default), averages over precession-timescale variations,
rms, compute the root-mean-square average, which is faster to evaluate.
The required inputs are either (theta1,theta2,deltaphi,a,e,q,chi1,chi2) or (deltachi,kappa,chieff,a,e,q,chi1,chi2). In the latter case, deltachi can be omitted and will be resampled from its PN-weighted distribution. Some inputs will be ignored for some of the chip definitions (for instance the heuristic chip does not depend on either r or deltaphi), but dummy values are still requested for code consistency. The keywords arguments are only used for the averaged formulation and are passed to precession_average.
- Parameters:
- theta1: float, optional (default: None)
Angle between orbital angular momentum and primary spin.
- theta2: float, optional (default: None)
Angle between orbital angular momentum and secondary spin.
- deltaphi: float, optional (default: None)
Angle between the projections of the two spins onto the orbital plane.
- deltachi: float, optional (default: None)
Weighted spin difference.
- kappa: float, optional (default: None)
Asymptotic angular momentum.
- a: float, optional (default: None)
Semi-major axis.
- e: float, optional (default: None)
Eccentricity: 0<=e<1
- chieff: float, optional (default: None)
Effective spin.
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- chi1: float, optional (default: None)
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float, optional (default: None)
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- which: string, optional (default: “averaged”)
Select function behavior.
- **kwargs: unpacked dictionary, optional
Additional keyword arguments.
- Returns:
- chip: float
Effective precessing spin.
Examples
chip = precession.eval_chip(theta1=theta1,theta2=theta2,deltaphi=deltaphi,a=a,e=e,q=q,chi1=chi1,chi2=chi2,which="heuristic")
chip = precession.eval_chip(theta1=theta1,theta2=theta2,deltaphi=deltaphi,a=a,e=e,q=q,chi1=chi1,chi2=chi2,which="generalized")
chip = precession.eval_chip(theta1=theta1,theta2=theta2,deltaphi=deltaphi,a=a,e=e,q=q,chi1=chi1,chi2=chi2,which="averaged")
chip = precession.eval_chip(theta1=theta1,theta2=theta2,deltaphi=deltaphi,a=a,e=e,q=q,chi1=chi1,chi2=chi2,which="rms")
chip = precession.eval_chip(deltachi=deltachi,kappa=kappa,a=a,e=e,chieff=chieff,q=q,chi1=chi1,chi2=chi2,which="heuristic")
chip = precession.eval_chip(deltachi=deltachi,kappa=kappa,a=a,e=e,chieff=chieff,q=q,chi1=chi1,chi2=chi2,which="generalized")
chip = precession.eval_chip(deltachi=deltachi,kappa=kappa,a=a,e=e,chieff=chieff,q=q,chi1=chi1,chi2=chi2,which="averaged")
chip = precession.eval_chip(deltachi=deltachi,kappa=kappa,a=a,e=e,chieff=chieff,q=q,chi1=chi1,chi2=chi2,which="rms")
- eccentricity.eval_chip_averaged(kappa, a, e, chieff, q, chi1, chi2, **kwargs)[source]
Averaged definition of the effective precessing spin, see arxiv:2011.11948. Additional keywords arguments are passed to precession_average.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- **kwargs: unpacked dictionary, optional
Additional keyword arguments.
- Returns:
- chip: float
Effective precessing spin.
Examples
chip = precession.eval_chip_averaged(kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.eval_chip_generalized(theta1, theta2, deltaphi, q, chi1, chi2)[source]
Generalized definition of the effective precessing spin chip, see arxiv:2011.11948. This definition retains all variations on the precession timescale.
- Parameters:
- theta1: float
Angle between orbital angular momentum and primary spin.
- theta2: float
Angle between orbital angular momentum and secondary spin.
- deltaphi: float
Angle between the projections of the two spins onto the orbital plane.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- chip: float
Effective precessing spin.
Examples
chip = precession.eval_chip_generalized(theta1,theta2,deltaphi,q,chi1,chi2)
- eccentricity.eval_chip_heuristic(theta1, theta2, q, chi1, chi2)[source]
Heuristic definition of the effective precessing spin chip (Schmidt et al 2015), see arxiv:2011.11948. This definition inconsistently averages over some, but not all, variations on the precession timescale.
- Parameters:
- theta1: float
Angle between orbital angular momentum and primary spin.
- theta2: float
Angle between orbital angular momentum and secondary spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- chip: float
Effective precessing spin.
Examples
chip = precession.eval_chip_heuristic(theta1,theta2,q,chi1,chi2)
- eccentricity.eval_chip_rms(kappa, a, e, chieff, q, chi1, chi2)[source]
Root-mean-square definition of the effective precessing spin. This is much faster to evaluate than the plain average.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- chip: float
Effective precessing spin.
Examples
chip = precession.eval_chip_rms(kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.eval_cosdeltaphi(deltachi, kappa, a, e, chieff, q, chi1, chi2)[source]
Cosine of the angle between the projections of the two spins onto the orbital plane.
- Parameters:
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- cosdeltaphi: float
Cosine of the angle between the projections of the two spins onto the orbital plane.
Examples
cosdeltaphi = precession.eval_cosdeltaphi(deltachi,kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.eval_costheta1(deltachi, chieff, q, chi1)[source]
Cosine of the angle between the orbital angular momentum and the spin of the primary black hole.
- Parameters:
- deltachi: float
Weighted spin difference.
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- Returns:
- costheta1: float
Cosine of the angle between orbital angular momentum and primary spin.
Examples
costheta1 = precession.eval_costheta1(deltachi,chieff,q,chi1)
- eccentricity.eval_costheta12(theta1=None, theta2=None, deltaphi=None, deltachi=None, kappa=None, chieff=None, q=None, chi1=None, chi2=None)[source]
Cosine of the angle between the two spins. Valid inputs are either (theta1,theta2,deltaphi) or (deltachi,kappa,chieff,q,chi1,chi2).
- Parameters:
- theta1: float, optional (default: None)
Angle between orbital angular momentum and primary spin.
- theta2: float, optional (default: None)
Angle between orbital angular momentum and secondary spin.
- deltaphi: float, optional (default: None)
Angle between the projections of the two spins onto the orbital plane.
- deltachi: float, optional (default: None)
Weighted spin difference.
- kappa: float, optional (default: None)
Asymptotic angular momentum.
- chieff: float, optional (default: None)
Effective spin.
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- chi1: float, optional (default: None)
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float, optional (default: None)
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- costheta12: float
Cosine of the angle between the two spins.
Examples
costheta12 = precession.eval_costheta12(theta1=theta1,theta2=theta2,deltaphi=deltaphi)
costheta12 = precession.eval_costheta12(deltachi=deltachi,kappa=kappa,chieff=chieff,q=q,chi1=chi1,chi2=chi2)
- eccentricity.eval_costheta2(deltachi, chieff, q, chi2)[source]
Cosine of the angle between the orbital angular momentum and the spin of the secondary black hole.
- Parameters:
- deltachi: float
Weighted spin difference.
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- costheta2: float
Cosine of the angle between orbital angular momentum and secondary spin.
Examples
costheta2 = precession.eval_costheta2(deltachi,chieff,q,chi2)
- eccentricity.eval_costhetaL(deltachi, kappa, a, e, chieff, q)[source]
Cosine of the angle betwen the orbital angular momentum and the total angular momentum.
- Parameters:
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- Returns:
- costhetaL: float
Cosine of the angle betwen orbital angular momentum and total angular momentum.
Examples
costhetaL = precession.eval_costhetaL(deltachi,kappa,a,e,chieff,q)
- eccentricity.eval_cyclesign(ddeltachidt=None, deltaphi=None, Lvec=None, S1vec=None, S2vec=None)[source]
- Evaluate if the input parameters are in the first of the second half of a precession cycle. We refer to this as the ‘sign’ of a precession cycle, defined as +1 if deltachi is increasing and -1 deltachi is decreasing. Valid inputs are one and not more of the following:
dSdt
deltaphi
Lvec, S1vec, S2vec.
- Parameters:
- ddeltachidt: float, optional (default: None)
Time derivative of the total spin.
- deltaphi: float, optional (default: None)
Angle between the projections of the two spins onto the orbital plane.
- Lvec: array, optional (default: None)
Cartesian vector of the orbital angular momentum.
- S1vec: array, optional (default: None)
Cartesian vector of the primary spin.
- S2vec: array, optional (default: None)
Cartesian vector of the secondary spin.
- Returns:
- cyclesign: integer
Sign (either +1 or -1) to cover the two halves of a precesion cycle.
Examples
cyclesign = precession.eval_cyclesign(ddeltachidt=ddeltachidt)
cyclesign = precession.eval_cyclesign(deltaphi=deltaphi)
cyclesign = precession.eval_cyclesign(Lvec=Lvec,S1vec=S1vec,S2vec=S2vec)
- eccentricity.eval_delta_omega(kappa, a, e, chieff, q, chi1, chi2, precomputedroots=None)[source]
Variation of the precession frequency of deltachi as it oscillates from deltachi- to deltachi+ back to deltachi-. This of the five phenomenological parameters from arxiv:2103.03894.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- delta_omega: float
Precession frequency variation due to nutation.
Examples
delta_omega = precession.eval_delta_omega(kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.eval_delta_theta(kappa, a, e, chieff, q, chi1, chi2, precomputedroots=None)[source]
Nutation amplitude of deltachi as it oscillates from deltachi- to deltachi+ back to deltachi-. This of the five phenomenological parameters from arxiv:2103.03894.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- delta_theta: float
Nutation amplitude.
Examples
delta_theta = precession.eval_delta_theta(kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.eval_deltachi(theta1, theta2, q, chi1, chi2)[source]
Weighted spin difference.
- Parameters:
- theta1: float
Angle between orbital angular momentum and primary spin.
- theta2: float
Angle between orbital angular momentum and secondary spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- deltachi: float
Weighted spin difference.
Examples
deltachi = precession.eval_deltachi(theta1,theta2,q,chi1,chi2)
- eccentricity.eval_deltachiinf(kappa, chieff, q, chi1, chi2)[source]
Large-separation limit of the weighted spin difference.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- deltachi: float
Weighted spin difference.
Examples
deltachi = precession.eval_deltachiinf(kappa,chieff,q,chi1,chi2)
- eccentricity.eval_deltaphi(deltachi, kappa, a, e, chieff, q, chi1, chi2, cyclesign=1)[source]
Angle between the projections of the two spins onto the orbital plane.
- Parameters:
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- cyclesign: integer, optional (default: 1)
Sign (either +1 or -1) to cover the two halves of a precesion cycle.
- Returns:
- deltaphi: float
Angle between the projections of the two spins onto the orbital plane.
Examples
deltaphi = precession.eval_deltaphi(deltachi,kappa,a,e,chieff,q,chi1,chi2,cyclesign=1)
- eccentricity.eval_e(L=None, u=None, uc=None, a=None, q=None)[source]
Orbital eccentricity of the binary. Valid inputs are either (L,a,q) or (u,uc,q).
- Parameters:
- L: float, optional (default: None)
Magnitude of the Newtonian orbital angular momentum.
- u: float, optional (default: None)
Compactified separation 1/(2L).
- uc: float, optional (default: None)
Circular compactified separation uc=u(e=0).
- a: float, optional (default: None)
Binary semi-major axis.
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- Returns:
- e: float
Binary eccentricity.
- eccentricity.eval_eta(q)[source]
Symmetric mass ratio eta = m1*m2/(m1+m2)^2 = q/(1+q)^2.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- Returns:
- eta: float
Symmetric mass ratio 0<=eta<=1/4.
Examples
eta = precession.eval_eta(q)
- eccentricity.eval_kappa(theta1=None, theta2=None, deltaphi=None, J=None, a=None, e=None, q=None, chi1=None, chi2=None)[source]
Asymptotic angular momentum. Provide either (theta1,theta2,deltaphi,a,e,q,chi1,chhi2) or (J,a,e,q).
- Parameters:
- theta1: float, optional (default: None)
Angle between orbital angular momentum and primary spin.
- theta2: float, optional (default: None)
Angle between orbital angular momentum and secondary spin.
- deltaphi: float, optional (default: None)
Angle between the projections of the two spins onto the orbital plane.
- J: float, optional (default: None)
Magnitude of the total angular momentum.
- a: float, optional (default: None)
Semi-major axis.
- e: float, optional (default: None)
Eccentricity: 0<=e<1
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- chi1: float, optional (default: None)
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float, optional (default: None)
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- kappa: float
Asymptotic angular momentum.
Examples
kappa = precession.eval_kappa(theta1=theta1,theta2=theta2,deltaphi=deltaphi,a=a,e=e,q=q,chi1=chi1,chi2=chi2)
kappa = precession.eval_kappa(J=J,a=a,e=e,q=q)
- eccentricity.eval_m1(q)[source]
Mass of the heavier black hole in units of the total mass.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- Returns:
- m1: float
Mass of the primary (heavier) black hole.
Examples
m1 = precession.eval_m1(q)
- eccentricity.eval_m2(q)[source]
Mass of the lighter black hole in units of the total mass.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- Returns:
- m2: float
Mass of the secondary (lighter) black hole.
Examples
m2 = precession.eval_m2(q)
- eccentricity.eval_nutation_freq(kappa, a, e, chieff, q, chi1, chi2, precomputedroots=None)[source]
Nutation frequency of deltachi as it oscillates from deltachi- to deltachi+ back to deltachi-. This of the five phenomenological parameters from arxiv:2103.03894.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- littleomega: float
Squared time derivative of the weighted spin difference.
Examples
littleomega = precession.eval_nutation_freq(kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.eval_phiL(deltachi, kappa, a, e, chieff, q, chi1, chi2, cyclesign=1, precomputedroots=None)[source]
Evaluate the azimuthal angle of L in a plane orthogonal to J along the precession cycle.
- Parameters:
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- cyclesign: integer, optional (default: 1)
Sign (either +1 or -1) to cover the two halves of a precesion cycle.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- phiL: float
Azimuthal angle spanned by L about J.
Examples
phiL = precession.eval_phiL(deltachi,kappa,a,e,chieff,q,chi1,chi2,cyclesign=1)
- eccentricity.eval_q(m1, m2)[source]
Mass ratio, 0 < q = m2/m1 < 1.
- Parameters:
- m1: float
Mass of the primary (heavier) black hole.
- m2: float
Mass of the secondary (lighter) black hole.
- Returns:
- q: float
Mass ratio: 0<=q<=1.
Examples
q = precession.eval_q(m1,m2)
- eccentricity.eval_tau(kappa, a, e, chieff, q, chi1, chi2, precomputedroots=None, return_psiperiod=False, donotnormalize=False)[source]
Evaluate the nutation period. Setting return_psiperiod=True return the phase-period instead of the time period (i.e. returns tau/2K(m) insteaf of tau). This is useful to avoid the evaluation of an elliptic integral when it’s not needed. Setting donotnormalize=True is equivalent to mathcalA=1.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- MISSING: COULD NOT BUILD, optional (default: False)
FILL MANUALLY.
- donotnormalize: boolean, optional (default: False)
If True omit the numerical prefactor.
- Returns:
- tau: float
Nutation period.
Examples
tau = precession.eval_tau(kappa,a,e,chieff,q,chi1,chi2)
tau = precession.eval_tau(kappa,a,e,chieff,q,chi1,chi2,return_psiperiod=True)
- eccentricity.eval_theta1(deltachi, chieff, q, chi1)[source]
Angle between the orbital angular momentum and the spin of the primary black hole.
- Parameters:
- deltachi: float
Weighted spin difference.
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- Returns:
- theta1: float
Angle between orbital angular momentum and primary spin.
Examples
theta1 = precession.eval_theta1(deltachi,chieff,q,chi1)
- eccentricity.eval_theta12(theta1=None, theta2=None, deltaphi=None, deltachi=None, kappa=None, chieff=None, q=None, chi1=None, chi2=None)[source]
Angle between the two spins. Valid inputs are either (theta1,theta2,deltaphi) or (deltachi,kappa,chieff,q,chi1,chi2).
- Parameters:
- theta1: float, optional (default: None)
Angle between orbital angular momentum and primary spin.
- theta2: float, optional (default: None)
Angle between orbital angular momentum and secondary spin.
- deltaphi: float, optional (default: None)
Angle between the projections of the two spins onto the orbital plane.
- deltachi: float, optional (default: None)
Weighted spin difference.
- kappa: float, optional (default: None)
Asymptotic angular momentum.
- chieff: float, optional (default: None)
Effective spin.
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- chi1: float, optional (default: None)
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float, optional (default: None)
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- theta12: float
Angle between the two spins.
Examples
theta12 = precession.eval_theta12(theta1=theta1,theta2=theta2,deltaphi=deltaphi)
theta12 = precession.eval_theta12(deltachi=deltachi,kappa=kappa,chieff=chieff,q=q,chi1=chi1,chi2=chi2)
- eccentricity.eval_theta2(deltachi, chieff, q, chi2)[source]
Angle between the orbital angular momentum and the spin of the secondary black hole.
- Parameters:
- deltachi: float
Weighted spin difference.
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- theta2: float
Angle between orbital angular momentum and secondary spin.
Examples
theta2 = precession.eval_theta2(deltachi,chieff,q,chi2)
- eccentricity.eval_thetaL(deltachi, kappa, a, e, chieff, q)[source]
Angle betwen the orbital angular momentum and the total angular momentum.
- Parameters:
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- Returns:
- thetaL: float
Angle betwen orbital angular momentum and total angular momentum.
Examples
thetaL = precession.eval_thetaL(deltachi,kappa,a,e,chieff,q)
- eccentricity.eval_u(a, e, q)[source]
Compactified orbital separation.
- Parameters:
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- q: float
Mass ratio: 0<=q<=1.
- Returns:
- u: float
Compactified separation 1/(2L).
Examples
u = precession.eval_u(r,q)
- eccentricity.eval_v(a, e)[source]
Newtonian orbital velocity of the binary.
- Parameters:
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- Returns:
- v: float
Newtonian orbital velocity.
Examples
v = precession.eval_v(r)
- eccentricity.implicit(u, uc)[source]
LHS of eq. (28) in Fumagalli & Gerosa, arXiv:2310.16893. This is used to compute the implicit function u(uc) in the inspiral_precav function. :param u: Compactified separation 1/(2L). :type u: float :param uc: Circular compactified separation 1/(2L(e=0)). :type uc: float
- eccentricity.inspiral(*args, which=None, **kwargs)[source]
- Generic wrapper for binary inspirals. The keyword which selects between the different approximation:
If which=precession use precession-averaged inspiral (i.e. run inspiral_precav)
If which=orbit use orbit-averaged inspiral (i.e. run inspiral_orbav)
If which=hybrid use hybrid inspiral (i.e. run inspiral_hybrid)
All other aguments and keyword arguments are passed on to the relevant inspiral function; see their docs for details.
- Parameters:
- **args: unpacked dictionary, optional
Additional arguments.
- which: string, optional (default: None)
Select function behavior.
- **kwargs: unpacked dictionary, optional
Additional keyword arguments.
- Returns:
- outputs: dictionary, optional
Set of outputs.
Examples
outputs = precession.inspiral(*args,which='precession',**kwargs)
outputs = precession.inspiral(*args,which='orbit',**kwargs)
outputs = precession.inspiral(*args,which='hybrid',**kwargs)
- eccentricity.inspiral_hybrid(theta1=None, theta2=None, deltaphi=None, deltachi=None, kappa=None, a=None, aswitch=None, e=None, uc=None, ucswitch=None, u=None, chieff=None, q=None, chi1=None, chi2=None, requested_outputs=None, **odeint_kwargs)[source]
Perform hybrid inspirals, i.e. evolve the binary at large separation with a pression-averaged evolution and at small separation with an orbit-averaged evolution, matching the two. The variables q, chi1, and chi2 must always be provided. The integration range must be specified using either a or uc (and not both); provide also ucswitch and aswitch consistently. Either a of uc needs to be arrays with lenght >=1, where e.g. a[0] corresponds to the initial condition and a[1:] corresponds to the location where outputs are returned. It does not work at past time infinity if e !=0. The function is vectorized: evolving N multiple binaries with M outputs requires kappainitial, chieff, q, chi1, chi2 to be of shape (N,) and u of shape (M,N). The initial conditions must be specified in terms of one an only one of the following:
theta1,theta2, and deltaphi (but note that deltaphi is not necessary if integrating from infinite separation).
kappa, chieff.
The desired outputs can be specified with a list e.g. requested_outputs=[‘theta1’,’theta2’,’deltaphi’]. All the available variables are returned by default. These are: [‘theta1’, ‘theta2’, ‘deltaphi’, ‘deltachi’, ‘kappa’, ‘a’,’e’, ‘uc’,’u’, ‘deltachiminus’, ‘deltachiplus’, ‘deltachi3’, ‘chieff’, ‘q’, ‘chi1’, ‘chi2’]. The flag enforce allows checking the consistency of the input variables. Additional keywords arguments are passed to scipy.integrate.odeint after some custom-made default settings.
- Parameters:
- theta1: float, optional (default: None)
Angle between orbital angular momentum and primary spin.
- theta2: float, optional (default: None)
Angle between orbital angular momentum and secondary spin.
- deltaphi: float, optional (default: None)
Angle between the projections of the two spins onto the orbital plane.
- deltachi: float, optional (default: None)
Weighted spin difference.
- kappa: float, optional (default: None)
Asymptotic angular momentum.
- a: float, optional (default: None)
Binary semi-major axis.
- e: float, optional (default: None)
Binary eccentricity: 0<=e<1 .
- aswitch: float, optional (default: None)
Matching separation between the precession- and orbit-averaged chunks.
- eswitch: float, optional (default: None)
Matching separation between the precession- and orbit-averaged chunks.
- u: float, optional (default: None)
Compactified separation 1/(2L).
- uc: float, optional (default: None)
Compactified circular separation 1/(2L).
- uswitch: float, optional (default: None)
Matching compactified separation between the precession- and orbit-averaged chunks.
- chieff: float, optional (default: None)
Effective spin.
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- chi1: float, optional (default: None)
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float, optional (default: None)
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- requested_outputs: list, optional (default: None)
Set of outputs.
- **odeint_kwargs: unpacked dictionary, optional
Additional keyword arguments.
- Returns:
- outputs: dictionary
Set of outputs.
Examples
outputs = precession.inspiral_hybrid(theta1=theta1,theta2=theta2,deltaphi=deltaphi,a=a,e=e,q=q,chi1=chi1,chi2=chi2)
outputs = precession.inspiral_hybrid(theta1=theta1,theta2=theta2,deltaphi=deltaphi,uc=uc,u=u,q=q,chi1=chi1,chi2=chi2)
outputs = precession.inspiral_hybrid(kappa,a=a,e=e,chieff=chieff,q=q,chi1=chi1,chi2=chi2)
outputs = precession.inspiral_hybrid(kappa,uc=uc,u=u,chieff=chieff,q=q,chi1=chi1,chi2=chi2)
- eccentricity.inspiral_orbav(theta1=None, theta2=None, deltaphi=None, Lh=None, S1h=None, S2h=None, delta_lambda=0, deltachi=None, kappa=None, a=None, e=None, uc=None, u=None, chieff=None, q=None, chi1=None, chi2=None, cyclesign=1, PNorderpre=[0, 0.5], PNorderrad=[0, 1, 1.5, 2, 2.5, 3, 3.5], requested_outputs=None, **odeint_kwargs)[source]
Perform precession-averaged inspirals. The variables q, chi1, and chi2 must always be provided. The integration range must be specified using either (a,e) or (uc,u) (and not both). These need to be arrays with lenght >=1, where e.g. a[0] corresponds to the initial condition and a[1:] corresponds to the location where outputs are returned. The function is vectorized: evolving N multiple binaries with M outputs requires kappainitial, chieff, q, chi1, chi2 to be of shape (N,) and uc of shape (M,N). The initial conditions must be specified in terms of one an only one of the following:
Lh, S1h, and S2h
theta1,theta2, and deltaphi.
deltachi, kappa, chieff, cyclesign.
The desired outputs can be specified with a list e.g. requested_outputs=[‘theta1’,’theta2’,’deltaphi’]. All the available variables are returned by default. These are: [‘theta1’, ‘theta2’, ‘deltaphi’, ‘deltachi’, ‘kappa’, ‘r’, ‘u’, ‘deltachiminus’, ‘deltachiplus’, ‘deltachi3’, ‘chieff’, ‘q’, ‘chi1’, ‘chi2’]. The flag enforce allows checking the consistency of the input variables. Additional keywords arguments are passed to scipy.integrate.odeint after some custom-made default settings.
- Parameters:
- theta1: float, optional (default: None)
Angle between orbital angular momentum and primary spin.
- theta2: float, optional (default: None)
Angle between orbital angular momentum and secondary spin.
- deltaphi: float, optional (default: None)
Angle between the projections of the two spins onto the orbital plane.
- Lh: array, optional (default: None)
Direction of the orbital angular momentum, unit vector.
- S1h: array, optional (default: None)
Direction of the primary spin, unit vector.
- S2h: array, optional (default: None)
Direction of the secondary spin, unit vector.
- deltachi: float, optional (default: None)
Weighted spin difference.
- kappa: float, optional (default: None)
Asymptotic angular momentum.
- a: float, optional (default: None)
Binary semi-major axis.
- e: float, optional (default: None)
Eccentricity 0<=e<1.
- uc: float, optional (default: None)
Circualr compactified separation 1/(2L(e=0)).
- u: float, optional (default: None)
Compactified separation 1/(2L).
- chieff: float, optional (default: None)
Effective spin.
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- chi1: float, optional (default: None)
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float, optional (default: None)
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- cyclesign: integer, optional (default: +1)
Sign (either +1 or -1) to cover the two halves of a precesion cycle.
- PNorderpre: array (default: [0,0.5])
PN orders considered in the spin-precession equations.
- PNorderrad: array (default: [0,0.5])
PN orders considered in the radiation-reaction equation.
- requested_outputs: list, optional (default: None)
Set of outputs.
- **odeint_kwargs: unpacked dictionary, optional
Additional keyword arguments.
- Returns:
- outputs: dictionary
Set of outputs.
Examples
outputs = precession.inspiral_orbav(Lh=Lh,S1h=S1h,S2h=S2h,a=a,e=e,q=q,chi1=chi1,chi2=chi2)
outputs = precession.inspiral_orbav(Lh=Lh,S1h=S1h,S2h=S2h,uc=uc,u=u,q=q,chi1=chi1,chi2=chi2)
outputs = precession.inspiral_orbav(theta1=theta1,theta2=theta2,deltaphi=deltaphi,a=a,e=e,q=q,chi1=chi1,chi2=chi2)
outputs = precession.inspiral_orbav(theta1=theta1,theta2=theta2,deltaphi=deltaphi,uc=uc,u=u,q=q,chi1=chi1,chi2=chi2)
outputs = precession.inspiral_orbav(deltachi=deltachi,kappa=kappa,a=a,e=e,chieff=chieff,q=q,chi1=chi1,chi2=chi2)
- eccentricity.inspiral_precav(theta1=None, theta2=None, deltaphi=None, deltachi=None, kappa=None, a=None, e=None, u=None, uc=None, chieff=None, q=None, chi1=None, chi2=None, requested_outputs=None, enforce=False, **odeint_kwargs)[source]
The integration range must be specified using either (a,e) or (uc, u), (and not both). These need to be arrays with lenght >=1, where e.g. a[0] corresponds to the initial condition and a[1:] corresponds to the location where outputs are returned. Do not go to past time infinity with e!=0. The function is vectorized: evolving N multiple binaries with M outputs requires kappainitial, chieff, q, chi1, chi2 to be of shape (N,) and u of shape (M,N). The initial conditions must be specified in terms of one an only one of the following:
theta1,theta2, and deltaphi.
kappa, chieff.
The desired outputs can be specified with a list e.g. requested_outputs=[‘theta1’,’theta2’,’deltaphi’]. All the available variables are returned by default. These are: [‘theta1’, ‘theta2’, ‘deltaphi’, ‘deltachi’, ‘kappa’, ‘a’, ‘e’, ‘u’, ‘deltachiminus’, ‘deltachiplus’, ‘deltachi3’, ‘chieff’, ‘q’, ‘chi1’, ‘chi2’]. The flag enforce allows checking the consistency of the input variables. Additional keywords arguments are passed to scipy.integrate.odeint after some custom-made default settings.
- Parameters:
- theta1: float, optional (default: None)
Angle between orbital angular momentum and primary spin.
- theta2: float, optional (default: None)
Angle between orbital angular momentum and secondary spin.
- deltaphi: float, optional (default: None)
Angle between the projections of the two spins onto the orbital plane.
- kappa: float, optional (default: None)
Asymptotic angular momentum.
- r: float, optional (default: None)
Binary separation.
- u: float, optional (default: None)
Compactified separation 1/(2L).
- chieff: float, optional (default: None)
Effective spin.
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- chi1: float, optional (default: None)
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float, optional (default: None)
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- requested_outputs: list, optional (default: None)
Set of outputs.
- enforce: boolean, optional (default: False)
If True raise errors, if False raise warnings.
- **odeint_kwargs: unpacked dictionary, optional
Additional keyword arguments.
- Returns:
- outputs: dictionary
Set of outputs.
Examples
outputs = precession.inspiral_precav(theta1=theta1,theta2=theta2,deltaphi=deltaphi,a=a,e=e,q=q,chi1=chi1,chi2=chi2)
outputs = precession.inspiral_precav(theta1=theta1,theta2=theta2,deltaphi=deltaphi,uc=uc,u=u,q=q,chi1=chi1,chi2=chi2)
outputs = precession.inspiral_precav(kappa,a=a,e=e,chieff=chieff,q=q,chi1=chi1,chi2=chi2)
outputs = precession.inspiral_precav(kappa,uc=uc,u=u,chieff=chieff,q=q,chi1=chi1,chi2=chi2)
- eccentricity.integrator_orbav(Lhinitial, S1hinitial, S2hinitial, delta_lambda, a, e, q, chi1, chi2, PNorderpre=[0, 0.5], PNorderrad=[0, 1, 1.5, 2, 2.5, 3], **odeint_kwargs)[source]
Integration of the systems of ODEs describing orbit-averaged inspirals. Additional keywords arguments are passed to scipy.integrate.odeint after some custom-made default settings.
- Parameters:
- Lhinitial: array
Initial direction of the orbital angular momentum, unit vector.
- S1hinitial: array
Initial direction of the primary spin, unit vector.
- S2hinitial: array
Initial direction of the secondary spin, unit vector.
- a: float
Bianary semi-major axis.
- efloat
Eccentricty: 0<=e<1. .
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- PNorderpre: array (default: [0,0.5])
PN orders considered in the spin-precession equations.
- PNorderrad: array (default: [0,0.5])
PN orders considered in the radiation-reaction equation.
- **odeint_kwargs: unpacked dictionary, optional
Additional keyword arguments.
- Returns:
- ODEsolution: array
Solution of the ODE.
Examples
ODEsolution = precession.integrator_orbavintegrator_orbav(Lhinitial,S1hinitial,S2hinitial,v,q,chi1,chi2)
- eccentricity.integrator_precav(kappainitial, u, chieff, q, chi1, chi2, **odeint_kwargs)[source]
Integration of ODE dkappa/du describing precession-averaged inspirals. Here u needs to be an array with lenght >=1, where u[0] corresponds to the initial condition and u[1:] corresponds to the location where outputs are returned. For past time infinity use u=0. The function is vectorized: evolving N multiple binaries with M outputs requires kappainitial, chieff, q, chi1, chi2 to be of shape (N,) and u of shape (M,N). Additional keywords arguments are passed to scipy.integrate.odeint after some custom-made default settings.
- Parameters:
- kappainitial: float
Initial value of the regularized momentum kappa.
- u: float
Compactified separation 1/(2L).
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- **odeint_kwargs: unpacked dictionary, optional
Additional keyword arguments.
- Returns:
- kappa: float
Asymptotic angular momentum.
Examples
kappa = precession.integrator_precav(kappainitial,u,chieff,q,chi1,chi2)
- eccentricity.intertial_ingredients(kappa, a, e, chieff, q, chi1, chi2)[source]
Numerical prefactors entering the precession frequency.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- bigC0: float
Prefactor in the OmegaL equation.
- bigCminus: float
Prefactor in the OmegaL equation.
- bigCplus: float
Prefactor in the OmegaL equation.
- bigRminus: float
Prefactor in the PhiL equation.
- bigRplus: float
Prefactor in the PhiL equation.
Examples
bigC0,bigCplus,bigCminus,bigRplus,bigRminus = precession.intertial_ingredients(kappa,a,e,chieff,q,chi1,chi2)
- eccentricity.inverseaffine(rescaled, low, up)[source]
Perform an inferse affine trasformation to scale a vector from some upper and lower limits: vec = low + rescaled*(up-low). If rescaled is defined in [0,1], the result will be defined in [low,up].
- Parameters:
- rescaled: float
Rescaled vector.
- low: float
Lower limit.
- up: float
Upper limit.
- Returns:
- vec: array
Vector in Cartesian coomponents.
Examples
vec = precession.inverseaffine(rescaled,low,up)
- eccentricity.ismonotonic(vec, which)[source]
- Check if an array is monotonic. The parameter which can takes the following values:
< check array is strictly increasing.
<= check array is increasing.
> check array is strictly decreasing.
>= check array is decreasing.
- Parameters:
- vec: array
Input array.
- which: string
Select function behavior.
- Returns:
- check: boolean
Result
Examples
check = ismonotonic(vec, which)
- eccentricity.isotropic_angles(N=1)[source]
Sample the angles theta1, theta2, and deltaphi from an isotropic distribution.
- Parameters:
- N: integer, optional (default: 1)
Number of samples.
- Returns:
- vec: array
Vector in Cartesian coomponents.
Examples
vec = precession.isotropic_angles(N=1)
- eccentricity.kappadiscriminant_coefficients(u, chieff, q, chi1, chi2)[source]
Coefficients of the quintic equation in kappa that defines the spin-orbit resonances.
- Parameters:
- u: float
Compactified separation 1/(2L).
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- coeff0: float
Coefficient to the x^0 term in polynomial.
- coeff1: float
Coefficient to the x^1 term in polynomial.
- coeff2: float
Coefficient to the x^2 term in polynomial.
- coeff3: float
Coefficient to the x^3 term in polynomial.
- coeff4: float
Coefficient to the x^4 term in polynomial.
- coeff5: float
Coefficient to the x^5 term in polynomial.
Examples
coeff5,coeff4,coeff3,coeff2,coeff1,coeff0 = precession.kappadiscriminant_coefficients(u,chieff,q,chi1,chi2)
- eccentricity.kappalimits(a=None, e=None, chieff=None, q=None, chi1=None, chi2=None, enforce=False, **kwargs)[source]
- Limits on the asymptotic angular momentum. The contraints considered depend on the inputs provided.
If r, q, chi1, chi2 are provided, returns the geometrical limits.
If r, chieff, q, chi1, and chi2 are provided, returns the spin-orbit resonances.
The boolean flag enforce raises an error in case the inputs are not compatible. Additional kwargs are passed to kapparesonances.
- Parameters:
- a: float, optional (default: None)
Semi-major axis.
- e: float, optional (default: None)
Eccentricity: 0<=e<1
- chieff: float, optional (default: None)
Effective spin.
- q: float, optional (default: None)
Mass ratio: 0<=q<=1.
- chi1: float, optional (default: None)
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float, optional (default: None)
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- enforce: boolean, optional (default: False)
Perform additional checks.
- **kwargs: unpacked dictionary, optional
Additional keyword arguments.
- Returns:
- kappamax: float
Maximum value of the asymptotic angular momentum kappa.
- kappamin: float
Minimum value of the asymptotic angular momentum kappa.
Examples
kappamin,kappamax = precession.kappalimits(a=a,e=e,q=q,chi1=chi1,chi2=chi2)
kappamin,kappamax = precession.kappalimits(a=a,e=e,chieff=chieff,q=q,chi1=chi1,chi2=chi2)
kappamin,kappamax = precession.kappalimits(a=a,e=e,chieff=chieff,q=q,chi1=chi1,chi2=chi2,enforce=True)
- eccentricity.kappalimits_geometrical(a, e, q, chi1, chi2)[source]
Limits in kappa conditioned on r, q, chi1, chi2.
- Parameters:
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- kappamax: float
Maximum value of the asymptotic angular momentum kappa.
- kappamin: float
Minimum value of the asymptotic angular momentum kappa.
Examples
kappamin,kappamax = precession.kappalimits_geometrical(r,q,chi1,chi2)
- eccentricity.kapparescaling(kappatilde, a, e, chieff, q, chi1, chi2)[source]
Compute kappa from the rescaled parameter 0<=kappatilde<=1.
- Parameters:
- kappatilde: float
Rescaled version of the asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- kappa: float
Asymptotic angular momentum.
Examples
kappa = precession.kapparescaling(kappatilde,a,e,chieff,q,chi1,chi2)
- eccentricity.kapparesonances(a, e, chieff, q, chi1, chi2, tol=0.0001)[source]
asymptotic angular momentum of the two spin-orbit resonances. The resonances minimizes and maximizes kappa for given values of r, chieff, q, chi1, chi2. The minimum corresponds to deltaphi=pi and the maximum corresponds to deltaphi=0.
- Parameters:
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- tol: float, optional (default: 1e-4)
Numerical tolerance, see source code for details.
- Returns:
- kappamax: float
Maximum value of the asymptotic angular momentum kappa.
- kappamin: float
Minimum value of the asymptotic angular momentum kappa.
Examples
kappamin,kappamax = precession.kapparesonances(u,chieff,q,chi1,chi2,tol=1e-4)
- eccentricity.morphology(kappa, a, e, chieff, q, chi1, chi2, simpler=False, precomputedroots=None)[source]
- Evaluate the spin morphology. Returns:
L0 for librating about deltaphi=0,
Lpi for librating about deltaphi=pi
C- for circulating from deltaphi=pi to deltaphi=0,
C+ for circulating from deltaphi=0 to deltaphi=pi.
If simpler=True, do not distinguish between the two circulating morphologies and return C for both.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- simpler: boolean, optional (default: False)
If True simplifies output.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- morph: string
Spin morphology.
Examples
morph = precession.morphology(kappa,a,e,chieff,q,chi1,chi2,simpler=False)
- eccentricity.norm_nested(x)[source]
Norm of 2D array of shape (N,3) along last axis.
- Parameters:
- xarray
Input array.
- Returns:
- narray
Norm of the input arrays.
Examples
n = norm_nested(x)
- eccentricity.normalize_nested(x)[source]
Normalize 2D array of shape (N,3) along last axis.
- Parameters:
- xarray
Input array.
- Returns:
- yarray
Normalized array.
Examples
y = normalize_nested(x)
- eccentricity.omegasq_aligned(a, e, q, chi1, chi2, which)[source]
Squared oscillation frequency of a given perturbed aligned-spin binary. The flag which needs to be set to uu for up-up, ud for up-down, du for down-up or dd for down-down where the term before (after) the hyphen refers to the spin of the heavier (lighter) black hole.
- Parameters:
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- which: string
Select function behavior.
- Returns:
- omegasq: float
Squared frequency.
Examples
omegasq = precession.omegasq_aligned(r,q,chi1,chi2,which)
- eccentricity.precession_average(kappa, a, e, chieff, q, chi1, chi2, func, *args, method='quadrature', Nsamples=10000.0)[source]
Average a generic function over a precession cycle. The function needs to have call: func(deltachi, *args). Keywords arguments are not supported. There are integration methods implemented:
method=’quadrature’ uses scipy.integrate.quad. This is set by default and should be preferred.
method=’montecarlo’ samples t(deltachi) and approximate the integral with a Monte Carlo sum. The number of samples can be specifed by Nsamples.
Additional keyword arguments are passed to scipy.integrate.quad.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- func: function
Function to precession-average.
- *args: tuple
Extra arguments to pass to func.
- method: string (default: ‘quadrature’)
Either ‘quadrature’ or ‘montecarlo’
- Nsamples: integer, optional (default: 1e4)
Number of Monte Carlo samples.
- Returns:
- func_av: float
Precession average of func.
Examples
func_av = precession.precession_average(kappa,a,e,chieff,q,chi1,chi2,func,*args,method='quadrature',Nsamples=1e4)
- eccentricity.rhs_orbav(allvars, a, q, m1, m2, eta, chi1, chi2, S1, S2, PNorderpre=[0, 0.5], PNorderrad=[0, 1, 1.5, 2, 2.5, 3])[source]
Right-hand side of the systems of ODEs describing orbit-averaged inspiral. The equations are reported in Sec 4A of Gerosa and Kesden, arXiv:1605.01067. The format is d[allvars]/dv=RHS where allvars=[Lhx,Lhy,Lhz,S1hx,S1hy,S1hz,S2hx,S2hy,S2hz,t], h indicates unit vectors, v is the orbital velocity, and t is time. This is an internal function used by the ODE integrator and is not array-compatible.
- Parameters:
- allvars: array
Packed ODE input variables.
- a: float
Newtonian orbital velocity.
- q: float
Mass ratio: 0<=q<=1.
- m1: float
Mass of the primary (heavier) black hole.
- m2: float
Mass of the secondary (lighter) black hole.
- eta: float
Symmetric mass ratio 0<=eta<=1/4.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- S1: float
Magnitude of the primary spin.
- S2: float
Magnitude of the secondary spin.
- PNorderpre: array (default: [0,0.5])
PN orders considered in the spin-precession equations.
- PNorderrad: array (default: [0,0.5])
PN orders considered in the radiation-reaction equation.
- Returns:
- RHS: float
Right-hand side.
Examples
RHS = precession.rhs_orbav(allvars,a,q,m1,m2,eta,chi1,chi2,S1,S2,PNorderpre=[0,0.5],PNorderrad=[0,1,1.5,2,2.5,3])
- eccentricity.rhs_precav(kappa, u, chieff, q, chi1, chi2)[source]
Right-hand side of the dkappa/du ODE describing precession-averaged inspiral. This is an internal function used by the ODE integrator and is not array-compatible.
- Parameters:
- kappa: float
Asymptotic angular momentum.
- u: float
Compactified separation 1/(2L).
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- RHS: float
Right-hand side.
Examples
RHS = precession.rhs_precav(kappa,u,chieff,q,chi1,chi2)
- eccentricity.roots_vec(p)[source]
Locate roots of polynomial using a vectorized version of numpy.roots. Credit to stackoverflow user ‘pv’ (https://stackoverflow.com/a/35853977) for some of this code, which we generalize to accomodate equations of different degrees.
- Parameters:
- p: float
Polynomial coefficients.
- Returns:
- roots: float
Polynomial roots.
Examples
roots = precession.roots_vec(p)
- eccentricity.rotate_nested(vec, align_zaxis, align_xzplane)[source]
Rotate a given vector vec to a frame such that the vector align_zaxis lies along z and the vector align_xzplane lies in the xz plane.
- Parameters:
- vecarray
Input array.
- align_zaxisarray
Reference array.
- align_xzplanearray
Reference array.
- Returns:
- vecrotarray
Rotated array.
Examples
vecrot = rotate_nested(vec, align_zaxis, align_xzplane)
- eccentricity.rupdown(q, chi1, chi2)[source]
The critical separations r_ud+/- marking the region of the up-down precessional instability.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- rudm: float
Inner orbital separation in the up-down instability.
- rudp: float
Outer orbital separation in the up-down instability.
Examples
rudp,rudm = precession.rupdown(q,chi1,chi2)
- eccentricity.sample_unitsphere(N=1)[source]
Sample points uniformly on a sphere of unit radius. Returns array of shape (N,3).
- Parameters:
- N: integer, optional (default: 1)
Number of samples.
- Returns:
- vec: array
Vector in Cartesian coomponents.
Examples
vec = sample_unitsphere(N=1)
- eccentricity.scalar_nested(k, x)[source]
Nested scalar product between a 1D and a 2D array.
- Parameters:
- kfloat
Input scalar.
- xarray
Input array.
- Returns:
- yarray
Scalar product array.
Examples
y = scalar_nested(k, x)
- eccentricity.tiler(thing, shaper)[source]
Repeat the quantity thing a numer of times give by the shape of shaper.
- Parameters:
- thing: float
Quantity to be repeated
- shaper: array
Quantity providing the shape.
- Returns:
- thing: float
Quantity to be repeated
Examples
thing = precession.tiler(thing,shaper)
- eccentricity.tofdeltachi(deltachi, kappa, a, e, chieff, q, chi1, chi2, cyclesign=1, precomputedroots=None)[source]
Time as a function of deltachi on the precessional timescale (without radiation reaction).
- Parameters:
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- chieff: float
Effective spin.
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- cyclesign: integer, optional (default: 1)
Sign (either +1 or -1) to cover the two halves of a precesion cycle.
- precomputedroots: array, optional (default: None)
Pre-computed output of deltachiroots for computational efficiency.
- Returns:
- t: float
Time.
Examples
t = precession.tofdeltachi(deltachi,kappa,a,e,chieff,q,chi1,chi2,cyclesign=1)
- eccentricity.updown_endpoint(q, chi1, chi2)[source]
Endpoint of the up-down precessional instability.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- deltaphi: float
Angle between the projections of the two spins onto the orbital plane.
- theta1: float
Angle between orbital angular momentum and primary spin.
- theta2: float
Angle between orbital angular momentum and secondary spin.
Examples
theta1,theta2,deltaphi = precession.updown_endpoint(q,chi1,chi2)
- eccentricity.vectors_to_Jframe(Lvec, S1vec, S2vec)[source]
Rotate vectors of the three momenta onto a frame where J is along z and L lies in the x-z plane.
- Parameters:
- Lvec: array
Cartesian vector of the orbital angular momentum.
- S1vec: array
Cartesian vector of the primary spin.
- S2vec: array
Cartesian vector of the secondary spin.
- Returns:
- Lvec: array
Cartesian vector of the orbital angular momentum.
- S1vec: array
Cartesian vector of the primary spin.
- S2vec: array
Cartesian vector of the secondary spin.
Examples
Lvec,S1vec,S2vec = precession.vectors_to_Jframe(Lvec,S1vec,S2vec)
- eccentricity.vectors_to_Lframe(Lvec, S1vec, S2vec)[source]
Rotate vectors of the three momenta onto a frame where L is along z and S1 lies in the x-z plane.
- Parameters:
- Lvec: array
Cartesian vector of the orbital angular momentum.
- S1vec: array
Cartesian vector of the primary spin.
- S2vec: array
Cartesian vector of the secondary spin.
- Returns:
- Lvec: array
Cartesian vector of the orbital angular momentum.
- S1vec: array
Cartesian vector of the primary spin.
- S2vec: array
Cartesian vector of the secondary spin.
Examples
Lvec,S1vec,S2vec = precession.vectors_to_Lframe(Lvec,S1vec,S2vec)
- eccentricity.vectors_to_angles(Lvec, S1vec, S2vec)[source]
Convert cartesian vectors (L,S1,S2) into angles (theta1,theta2,deltaphi).
- Parameters:
- Lvec: array
Cartesian vector of the orbital angular momentum.
- S1vec: array
Cartesian vector of the primary spin.
- S2vec: array
Cartesian vector of the secondary spin.
- Returns:
- deltaphi: float
Angle between the projections of the two spins onto the orbital plane.
- theta1: float
Angle between orbital angular momentum and primary spin.
- theta2: float
Angle between orbital angular momentum and secondary spin.
Examples
theta1,theta2,deltaphi = precession.vectors_to_angles(Lvec,S1vec,S2vec)
- eccentricity.vectors_to_conserved(Lvec, S1vec, S2vec, a, e, q, full_output=False)[source]
Convert vectors (L,S1,S2) to conserved quanties (deltachi,kappa,chieff).
- Parameters:
- Lvec: array
Cartesian vector of the orbital angular momentum.
- S1vec: array
Cartesian vector of the primary spin.
- S2vec: array
Cartesian vector of the secondary spin.
- q: float
Mass ratio: 0<=q<=1.
- full_output: boolean, optional (default: False)
Return additional outputs.
- Returns:
- chieff: float
Effective spin.
- cyclesign: integer, optional
Sign (either +1 or -1) to cover the two halves of a precesion cycle.
- deltachi: float
Weighted spin difference.
- kappa: float
Asymptotic angular momentum.
Examples
deltachi,kappa,chieff = precession.vectors_to_conserved(Lvec,S1vec,S2vec,q)
deltachi,kappa,chieff,cyclesign = precession.vectors_to_conserved(Lvec,S1vec,S2vec,q,full_output=True)
- eccentricity.widenutation_condition(a, e, q, chi1, chi2)[source]
- Conditions for wide nutation to take place. The returned flag which is
wide1 if wide nutation is allowed for the primary BH.
wide2 if wide nutation is allowed for the secondary BH.
nowide if wide nutation is not allowed.
- Parameters:
- a: float
Semi-major axis.
- e: float
Eccentricity: 0<=e<1
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- chieff: float
Effective spin.
- kappa: float
Asymptotic angular momentum.
- which: string
Select function behavior.
Examples
which,kappa,chieff = precession.widenutation_condition(r,q,chi1,chi2)
- eccentricity.widenutation_separation(q, chi1, chi2)[source]
The critical separation r_wide below which the binary component with smaller dimensionless spin may undergo wide nutations.
- Parameters:
- q: float
Mass ratio: 0<=q<=1.
- chi1: float
Dimensionless spin of the primary (heavier) black hole: 0<=chi1<=1.
- chi2: float
Dimensionless spin of the secondary (lighter) black hole: 0<=chi2<=1.
- Returns:
- r_wide: float
Orbital separation where wide nutations becomes possible.
Examples
r_wide = precession.widenutation_separation(q,chi1,chi2)
- eccentricity.wraproots(coefficientfunction, *args, **kwargs)[source]
Find roots of a polynomial given coefficients, ordered according to their real part. Complex roots are masked with nans. This is essentially a wrapper of numpy.roots.
- Parameters:
- coefficientfunction: callable
Function returning the polynomial coefficients ordered from highest to lowest degree.
- *args, **kwargs:
Parameters of coefficientfunction.
- Returns:
- sols: array
Roots of the polynomial.
Examples
sols = precession.wraproots(coefficientfunction, *args, **kwargs)