pint.models.binary_ddk.BinaryDDK

class pint.models.binary_ddk.BinaryDDK[source]

Bases: BinaryDD

Damour and Deruelle model with kinematics.

This extends the pint.models.binary_dd.BinaryDD model with “Shklovskii” and “Kopeikin” terms that account for the finite distance of the system from Earth, the finite size of the system, and the interaction of these with the proper motion.

From Kopeikin (1995) this includes \(\Delta_{\pi M}\) (Equation 17), the mixed annual-orbital parallax term, which changes \(a_1\) and \(\omega\) (delta_a1_parallax() and delta_omega_parallax()).

It does not include \(\Delta_{\pi P}\), the pure pulsar orbital parallax term (Equation 14).

From Kopeikin (1996) this includes apparent changes in \(\omega\), \(a_1\), and \(i\) due to the proper motion (delta_omega_proper_motion(), delta_a1_proper_motion(), delta_kin_proper_motion()) (Equations 8, 9, 10).

The actual calculations for this are done in pint.models.stand_alone_psr_binaries.DDK_model.DDKmodel.

It supports all the parameters defined in pint.models.pulsar_binary.PulsarBinary and pint.models.binary_dd.BinaryDD plus:

KIN

the inclination angle: \(i\)

KOM

the longitude of the ascending node, Kopeikin (1995) Eq 9: \(\Omega\)

K96

flag for Kopeikin binary model proper motion correction

It also removes:

SINI

use KIN instead

Note

This model defines KOM with reference to east, either equatorial or ecliptic depending on how the model is defined. KOM and KIN are defined in the Damour & Taylor (1992) convention (DT92), where:

KIN = 180 deg means the orbital angular momentum vector points toward the Earth, and KIN = 0 means the orbital angular momentum vector points away from the Earth.

KOM is 0 toward the East and increases clockwise on the sky; it is measured “East through North.”

Parameters supported:

Name / Aliases

Description

Kind

PB

Orbital period

d

PBDOT

Orbital period derivative respect to time

number

A1

Projected semi-major axis of pulsar orbit, ap*sin(i)

ls

A1DOT / XDOT

Derivative of projected semi-major axis, d[ap*sin(i)]/dt

ls / s

ECC / E

Eccentricity

number

EDOT

Eccentricity derivative respect to time

1 / s

T0

Epoch of periastron passage

d

OM

Longitude of periastron

deg

OMDOT

Rate of advance of periastron

deg / yr

M2

Companion mass

solMass

FB{number}

0th time derivative of frequency of orbit

1 / s

A0

DD model aberration parameter A0

s

B0

DD model aberration parameter B0

s

GAMMA

Time dilation & gravitational redshift

s

DR

Relativistic deformation of the orbit

number

DTH / DTHETA

Relativistic deformation of the orbit

number

KIN

Inclination angle

deg

KOM

The longitude of the ascending node

deg

K96

Flag for Kopeikin binary model proper motion correction

boolean

KINIAU

Inclination angle in the IAU convention

deg

KOMIAU

The longitude of the ascending node in the IAU convention

deg

SINI

Sine of inclination angle

number

References

  • Kopeikin (1995), ApJ, 439, L5 [1]

  • Kopeikin (1996), ApJ, 467, L93 [2]

  • Damour & Taylor (1992), Phys Rev D, 45, 1840 [3]

Methods

FBX_description(n)

FBX_unit(n)

add_param(param[, deriv_func, setup])

Add a parameter to the Component.

alternative_solutions()

Alternative Kopeikin solutions (potential local minima)

apply_units()

Apply units to parameter value.

binarymodel_delay(toas[, acc_delay])

Return the binary model independent delay call.

change_binary_epoch(new_epoch)

Change the epoch for this binary model.

check_required_params(required_params)

d_binary_delay_d_xxxx(toas, param, acc_delay)

Return the binary model delay derivatives.

get_params_of_type(param_type)

Get all the parameters in timing model for one specific type.

get_prefix_mapping_component(prefix)

Get the index mapping for the prefix parameters.

is_in_parfile(para_dict)

Check if this subclass included in parfile.

match_param_aliases(alias)

Return the parameter corresponding to this alias.

param_help()

Print help lines for all available parameters in model.

pb([t])

Return binary period and uncertainty (optionally evaluated at different times) regardless of binary model

print_par([format])

param format:

Parfile output format. PINT outputs the 'tempo', 'tempo2' and 'pint'

register_deriv_funcs(func, param)

Register the derivative function in to the deriv_func dictionaries.

remove_param(param)

Remove a parameter from the Component.

set_special_params(spcl_params)

setup()

Finalize construction loaded values.

update_binary_object(toas[, acc_delay])

Update stand alone binary's parameters and toas from PINT-facing object.

validate()

Validate parameters.

validate_toas(toas)

Check that this model component has TOAs where needed.

Attributes

PMLAT_DDK

Proper motion in latitude (Dec or ecliptic latitude)

PMLONG_DDK

Proper motion in longitude (RA or ecliptic longitude)

aliases_map

Return all the aliases and map to the PINT parameter name.

category

component_types

free_params_component

Return the free parameters in the component.

param_prefixs

register

property PMLONG_DDK

Proper motion in longitude (RA or ecliptic longitude)

property PMLAT_DDK

Proper motion in latitude (Dec or ecliptic latitude)

validate()[source]

Validate parameters.

alternative_solutions()[source]

Alternative Kopeikin solutions (potential local minima)

There are 4 potential local minima for a DDK model where a1dot is the same These are given by where Eqn. 8 in Kopeikin (1996) is equal to the best-fit value.

We first define the symmetry point where a1dot is zero (in equatorial coordinates):

\(KOM_0 = \tan^{-1} (\mu_{\delta} / \mu_{\alpha})\)

The solutions are then:

\((KIN, KOM)\)

\((KIN, 2KOM_0 - KOM - 180^{\circ})\)

\((180^{\circ}-KIN, KOM+180^{\circ})\)

\((180^{\circ}-KIN, 2KOM_0 - KOM)\)

All values will be between 0 and \(360^{\circ}\).

Returns:

tuple of (KIN,KOM) pairs for the four potential solutions

Return type:

tuple

add_param(param, deriv_func=None, setup=False)

Add a parameter to the Component.

The parameter is stored in an attribute on the Component object. Its name is also recorded in a list, self.params.

Parameters:
  • param (pint.models.Parameter) – The parameter to be added.

  • deriv_func (function) – Derivative function for parameter.

property aliases_map

Return all the aliases and map to the PINT parameter name.

This property returns a dictionary from the current in timing model parameters’ aliase to the pint defined parameter names. For the aliases of a prefixed parameter, the aliase with an existing prefix index maps to the PINT defined parameter name with the same index. Behind the scenes, the indexed parameter adds the indexed aliase to its aliase list.

apply_units()

Apply units to parameter value.

binarymodel_delay(toas, acc_delay=None)

Return the binary model independent delay call.

change_binary_epoch(new_epoch)

Change the epoch for this binary model.

T0 will be changed to the periapsis time closest to the supplied epoch, and the argument of periapsis (OM), eccentricity (ECC), and projected semi-major axis (A1 or X) will be updated according to the specified OMDOT, EDOT, and A1DOT or XDOT, if present.

Note that derivatives of binary orbital frequency higher than the first (FB2, FB3, etc.) are ignored in computing the new T0, even if present in the model. If high-precision results are necessary, especially for models containing higher derivatives of orbital frequency, consider re-fitting the model to a set of TOAs. The use of pint.simulation.make_fake_toas() and the pint.fitter.Fitter option track_mode="use_pulse_number" can make this extremely simple.

Parameters:

new_epoch (float MJD (in TDB) or astropy.Time object) – The new epoch value.

d_binary_delay_d_xxxx(toas, param, acc_delay)

Return the binary model delay derivatives.

property free_params_component

Return the free parameters in the component.

This function collects the non-frozen parameters.

Return type:

A list of free parameters.

get_params_of_type(param_type)

Get all the parameters in timing model for one specific type.

get_prefix_mapping_component(prefix)

Get the index mapping for the prefix parameters.

Parameters:

prefix (str) – Name of prefix.

Returns:

A dictionary with prefix parameter real index as key and parameter name as value.

Return type:

dict

is_in_parfile(para_dict)

Check if this subclass included in parfile.

Parameters:

para_dict (dictionary) – A dictionary contain all the parameters with values in string from one parfile

Returns:

Whether the subclass is included in the parfile.

Return type:

bool

match_param_aliases(alias)

Return the parameter corresponding to this alias.

Parameters:

alias (str) – Alias name.

Note

This function only searches the parameter aliases within the current component. If one wants to search the aliases in the scope of TimingModel, please use TimingModel.match_param_aliase().

param_help()

Print help lines for all available parameters in model.

pb(t=None)

Return binary period and uncertainty (optionally evaluated at different times) regardless of binary model

Parameters:

t (astropy.time.Time, astropy.units.Quantity, numpy.ndarray, float, int, str, optional) – Time(s) to evaluate period

Returns:

  • astropy.units.Quantity – Binary period

  • astropy.units.Quantity – Binary period uncertainty

print_par(format='pint')
Parameters:

format (str, optional) – Parfile output format. PINT outputs the ‘tempo’, ‘tempo2’ and ‘pint’ format. The defaul format is pint. Actual formatting done elsewhere.

Returns:

str

Return type:

formatted line for par file

register_deriv_funcs(func, param)

Register the derivative function in to the deriv_func dictionaries.

Parameters:
  • func (callable) – Calculates the derivative

  • param (str) – Name of parameter the derivative is with respect to

remove_param(param)

Remove a parameter from the Component.

Parameters:

param (str or pint.models.Parameter) – The parameter to remove.

setup()

Finalize construction loaded values.

update_binary_object(toas, acc_delay=None)

Update stand alone binary’s parameters and toas from PINT-facing object.

This function passes the PINT-facing object’s parameter values and TOAs to the stand-alone binary object. If the TOAs are not provided, it only updates the parameters not the TOAs.

Parameters:
  • toas (pint.toa.TOAs) – The TOAs that need to pass to the stand alone model.Default value is None. If toas is None, this function only updates the parameter value. If ‘acc_delay’ is not provided, the stand alone binary receives the standard barycentered TOAs.

  • acc_delay (numpy.ndarray) – If provided, TOAs will be corrected by provided acc_delay instead of the standard barycentering. The stand alone binary receives the input TOAs - acc_delay.

Notes

The values for obs_pos (the observatory position wrt the Solar System Barycenter) and psr_pos (the pulsar position wrt the Solar System Barycenter) are both computed in the same reference frame, ICRS or ECL depending on the model.

Warns:
  • If passing ‘None’ to ‘toa’ argument, the stand alone binary model will use

  • the TOAs were passed to it from last iteration (i.e. last barycentered

  • TOAs) or no TOAs for stand alone binary model at all. This behavior will

  • cause incorrect answers. Allowing the passing None to ‘toa’ argument is

  • for some lower level functions and tests. We do not recommend PINT

  • user to use it.

validate_toas(toas)

Check that this model component has TOAs where needed.