pint.models.stand_alone_psr_binaries.binary_generic.PSR_BINARY

class pint.models.stand_alone_psr_binaries.binary_generic.PSR_BINARY[source]

Bases: object

A base (generic) object for psr binary models.

In this class, a set of generally used binary parameters and several commonly used calculations are defined. For each binary model, the specific parameters and calculations are defined in the subclass.

A binary model takes the solar system barycentric time (ssb time) as input. When a binary model is instantiated, the parameters are set to the default values and input time is not initialized. The update those values, method .update_input() should be used.

Example of build a specific binary model class:

>>> from pint.models.stand_alone_psr_binaries.pulsar_binary import PSR_BINARY
>>> import numpy
>>> class foo(PSR_BINARY):
        def __init__(self):
            # This is to initialize the mother class attributes.
            super().__init__()
            self.binary_name = 'foo'
            # Add parameter that specific for my_binary, with default value and units
            self.param_default_value.update({'A0':0*u.second,'B0':0*u.second,
                                   'DR':0*u.Unit(''),'DTH':0*u.Unit(''),
                                   'GAMMA':0*u.second,})
            self.set_param_values() # This is to set all the parameters to attributes
            self.binary_delay_funcs += [self.foo_delay]
            self.d_binarydelay_d_par_funcs += [self.d_foo_delay_d_par]
            # If you have intermediate value in the calculation
            self.dd_interVars = ['er','eTheta','beta','alpha','Dre','Drep','Drepp',
                                 'nhat', 'TM2']
            self.add_inter_vars(self.dd_interVars)

        def foo_delay(self):
            pass

        def d_foo_delay_d_par(self):
            pass
>>> # To build a model instance
>>> binary_foo = foo()
>>> # binary_foo class has the default parameter value without toa input.
>>> # Update the toa input and parameters
>>> t = numpy.linspace(54200.0,55000.0,800)
>>> parameters_dict = {'A0':0.5,'ECC':0.01}
>>> binary_foo.update_input(t, parameters_dict)
>>> # Now the binary delay and derivatives can be computed.

To access the binary model class from pint platform, a pint pulsar binary wrapper is needed. See docstrings in the source code of pint/models/pulsar_binary class PulsarBinary.

Included general parameters: @param PB: Binary period [days] @param ECC: Eccentricity @param A1: Projected semi-major axis (lt-sec) @param A1DOT: Time-derivative of A1 (lt-sec/sec) @param T0: Time of periastron passage (barycentric MJD) @param OM: Omega (longitude of periastron) [deg] @param EDOT: Time-derivative of ECC [0.0] @param PBDOT: Time-derivative of PB [0.0] @param XPBDOT: Rate of change of orbital period minus GR prediction @param OMDOT: Time-derivative of OMEGA [0.0]

Intermediate variables calculation method are given here: Eccentric Anomaly E (not parameter ECC) Mean Anomaly M True Anomaly nu Eccentric ecc Longitude of periastron omega projected semi-major axis of orbit a1 TM2

Methods

Doppler()

E()

Eccentric Anomaly

M()

Orbit phase.

P()

Pobs_designmatrix(params)

TM2()

a1()

add_binary_params(parameter, defaultValue[, ...])

Add one parameter to the binary class.

add_inter_vars(interVars)

binary_delay()

Returns total pulsar binary delay.

compute_eccentric_anomaly(eccentricity, ...)

Solve the Kepler Equation, E - e * sin(E) = M

d_E_d_ECC()

d_E_d_EDOT()

d_E_d_T0()

Analytic derivative

d_E_d_par(par)

derivative for E respect to binary parameter.

d_M_d_par(par)

derivative for M respect to binary parameter.

d_Pobs_d_A1()

d_Pobs_d_A1DOT()

d_Pobs_d_ECC()

d_Pobs_d_EDOT()

d_Pobs_d_OM()

d_Pobs_d_OMDOT()

d_Pobs_d_P0()

d_Pobs_d_P1()

d_Pobs_d_PB()

d_Pobs_d_PBDOT()

d_Pobs_d_T0()

d_TM2_d_M2()

d_a1_d_A1()

d_a1_d_A1DOT()

d_a1_d_T0()

d_a1_d_par(par)

d_binarydelay_d_par(par)

Get the binary delay derivatives respect to parameters.

d_ecc_d_ECC()

d_ecc_d_EDOT()

d_ecc_d_T0()

d_nu_d_E()

d_nu_d_ECC()

Analytic calculation.

d_nu_d_EDOT()

d_nu_d_T0()

Analytic calculation.

d_nu_d_ecc()

d_nu_d_par(par)

derivative for nu respect to binary parameter.

d_omega_d_OM()

dOmega/dOM = 1

d_omega_d_OMDOT()

dOmega/dOMDOT = tt0

d_omega_d_T0()

dOmega/dPB = dnu/dPB*k+dk/dPB*nu

d_omega_d_par(par)

derivative for omega respect to user input Parameter.

d_pb_d_par(par)

derivative for pbprime respect to binary parameter.

delay_designmatrix(params)

ecc()

Calculate eccentricity with EDOT

get_tt0(barycentricTOA)

tt0 = barycentricTOA - T0

nu()

True anomaly (Ae)

omega()

orbits()

pb()

pbdot()

pbprime()

prtl_der(y, x)

Find the partial derivatives in binary model pdy/pdx

search_alias(parname)

set_param_values([valDict])

Set the parameters and assign values,

t0()

update_input(**updates)

Update the toas and parameters.

Attributes

T0

t

tt0

update_input(**updates)[source]

Update the toas and parameters.

set_param_values(valDict=None)[source]

Set the parameters and assign values,

If the valDict is not provided, it will set parameter as default value

add_binary_params(parameter, defaultValue, unit=False)[source]

Add one parameter to the binary class.

binary_delay()[source]

Returns total pulsar binary delay.

Returns:

Pulsar binary delay in the units of second

Return type:

float

d_binarydelay_d_par(par)[source]

Get the binary delay derivatives respect to parameters.

Parameters:

par (str) – Parameter name.

prtl_der(y, x)[source]

Find the partial derivatives in binary model pdy/pdx

Parameters:
  • y (str) – Name of variable to be differentiated

  • x (str) – Name of variable the derivative respect to

Returns:

The derivatives pdy/pdx

Return type:

np.array

compute_eccentric_anomaly(eccentricity, mean_anomaly)[source]

Solve the Kepler Equation, E - e * sin(E) = M

Parameters:
  • eccentricity (array_like) – Eccentricity of binary system

  • mean_anomaly (array_like) – Mean anomaly of the binary system

Returns:

The eccentric anomaly in radians, given a set of mean_anomalies in radians.

Return type:

array_like

get_tt0(barycentricTOA)[source]

tt0 = barycentricTOA - T0

ecc()[source]

Calculate eccentricity with EDOT

d_pb_d_par(par)[source]

derivative for pbprime respect to binary parameter.

Parameters:

par (string) – parameter name

Return type:

Derivative of M respect to par

M()[source]

Orbit phase.

d_M_d_par(par)[source]

derivative for M respect to binary parameter.

Parameters:

par (string) – parameter name

Return type:

Derivative of M respect to par

E()[source]

Eccentric Anomaly

d_E_d_T0()[source]

Analytic derivative

d(E-e*sinE)/dT0 = dM/dT0 dE/dT0(1-cosE*e)-de/dT0*sinE = dM/dT0 dE/dT0(1-cosE*e)+eDot*sinE = dM/dT0

d_E_d_par(par)[source]

derivative for E respect to binary parameter.

Parameters:

par (string) – parameter name

Return type:

Derivative of E respect to par

nu()[source]

True anomaly (Ae)

d_nu_d_T0()[source]

Analytic calculation.

dnu/dT0 = dnu/de*de/dT0+dnu/dE*dE/dT0 de/dT0 = -EDOT

d_nu_d_ECC()[source]

Analytic calculation.

dnu(e,E)/dECC = dnu/de*de/dECC+dnu/dE*dE/dECC de/dECC = 1 dnu/dPBDOT = dnu/dE*dE/dECC+dnu/de

d_nu_d_par(par)[source]

derivative for nu respect to binary parameter.

Parameters:

par (string) – parameter name

Return type:

Derivative of nu respect to par

d_omega_d_par(par)[source]

derivative for omega respect to user input Parameter.

Parameters:

par (string) – parameter name

Return type:

Derivative of omega respect to par

d_omega_d_OM()[source]

dOmega/dOM = 1

d_omega_d_OMDOT()[source]

dOmega/dOMDOT = tt0

d_omega_d_T0()[source]

dOmega/dPB = dnu/dPB*k+dk/dPB*nu