pint.models.stand_alone_psr_binaries.BT_piecewise.BTpiecewise
- class pint.models.stand_alone_psr_binaries.BT_piecewise.BTpiecewise(axis_store_initial=None, t=None, input_params=None)[source]
Bases:
BTmodel
This is a class independent from the PINT platform for pulsar BT piecewise binary model. It is a subclass of BTmodel which is a subclass of PSR_BINARY class defined in file binary_generic.py in the same directory. This class is designed for use with the PINT platform but can be used as an independent module for piecewise binary delay calculation. To interact with the PINT platform, a pulsar_binary wrapper is needed. See the source file pint/models/binary_piecewise.py.
Reference
The ‘BT’ binary model for the pulse period. Model as in: W.M. Smart, (1962), “Spherical Astronomy”, p35 Blandford & Teukolsky (1976), ApJ, 205, 580-591 :returns: * A piecewise bt binary model class with parameters, delay calculations and derivatives.
Example Session
———
>>import astropy.units as u
>>import numpy as np
>>binary_model=BTpiecewise()
>>param_dict = {‘T0’ (50000, ‘ECC’: 0.2})
>>binary_model.update_input(**param_dict)
>>t=np.linspace(50001.,60000.,10)*u.d
Adding binary parameters and piece ranges
>>binary_model.add_binary_params(‘T0X_0000’, 60000*u.d)
>>binary_model.add_binary_params(‘XR1_0000’, 50000*u.d)
>>binary_model.add_binary_params(‘XR2_0000’, 55000*u.d)
Can add more pieces here…
Overide default values values if desired
>>updates = {‘T0X_0000’ (60000.*u.d,’XR1_0000’:50000.*u.d,’XR2_0000’: 55000*u.d})
update the model with the piecewise parameter value(s) and piece ranges
>>binary_model.update_input(**updates)
Using pint’s get_model and loading this as a timing model and following the method described in ../binary_piecewise.py
sets _t multiple times during pint’s residual calculation
for simplicity we’re just going to set _t directly though this is not recommended.
>>setattr(binary_model,’_t’ ,t)
#here we call get_tt0 to get the “loaded toas” to interact with the pieces passed to the model earlier
#sets the attribute “T0X_per_toa” and/or “A1X_per_toa”, contains the piecewise parameter value that will be referenced
#for each toa future calculations
>>binary_model.get_tt0(t)
#For a piecewise T0, tt0 becomes a piecewise quantity, otherwise it is how it functions in BT_model.py.
#get_tt0 sets the attribute “T0X_per_toa” and/or “A1X_per_toa”.
#contains the piecewise parameter value that will be referenced for each toa future calculations
>>binary_model.T0X_per_toa
Information about any group can be found with the following
>>binary_model.piecewise_parameter_information
Order ([[Group index, Piecewise T0, Piecewise A1, Piece lower bound, Piece upper bound]])
Making sure a binary_model.tt0 exists
>>binary_model._tt0 = binary_model.get_tt0(binary_model._t)
Obtain piecewise BTdelay()
>>binary_model.BTdelay()
Methods
BTdelay
()Full BT model delay
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)Returns total pulsar binary delay.
compute_eccentric_anomaly
(eccentricity, ...)Solve the Kepler Equation, E - e * sin(E) = M
d_BTdelay_d_par
(par)d_E_d_ECC
()d_E_d_EDOT
()d_E_d_T0
()Analytic derivative
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)derivative for E respect to binary parameter.
d_M_d_par
(par)derivative for M respect to bianry 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_delayL1_d_A1
()d_delayL1_d_A1DOT
()d_delayL1_d_A1X
()d_delayL1_d_E
()d_delayL1_d_ECC
()d_delayL1_d_EDOT
()d_delayL1_d_GAMMA
()d_delayL1_d_OM
()d_delayL1_d_OMDOT
()d_delayL1_d_T0
()d_delayL1_d_T0X
()d_delayL1_d_par
(par)d_delayL2_d_A1
()d_delayL2_d_A1DOT
()d_delayL2_d_A1X
()d_delayL2_d_E
()d_delayL2_d_ECC
()d_delayL2_d_EDOT
()d_delayL2_d_GAMMA
()d_delayL2_d_OM
()d_delayL2_d_OMDOT
()d_delayL2_d_T0
()d_delayL2_d_T0X
()d_delayL2_d_par
(par)d_ecc_d_ECC
()d_ecc_d_EDOT
()d_ecc_d_T0
()d_nu_d_E
()Analytic calculation.
d_nu_d_EDOT
()Analytic calculation.
d_nu_d_ecc
()d_nu_d_par
(par)derivative for nu respect to binary parameter.
dOmega/dOM = 1
dOmega/dOMDOT = tt0
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.
delayL1
()First term of Blandford & Teukolsky (1976), ApJ, 205, 580-591, eq 2.33/ First left-hand term of W.M.
delayL2
()Second term of Blandford & Teukolsky (1976), ApJ, 205, 580-591, eq 2.33/ / Second left-hand term of W.M.
delayR
()Third term of Blandford & Teukolsky (1976), ApJ, 205, 580-591, eq 2.33 / Right-hand term of W.M.
delay_designmatrix
(params)ecc
()Calculate eccentricity with EDOT
Get the piecewise group boundaries from the dictionary of piecewise parameter information. :returns: list (length: number of pieces). Contains all pieces' lower edge. list (length: number of pieces). Contains all pieces' upper edge. :rtype: list.
get_tt0
(barycentricTOA)Finds (barycentricTOA - T0_x). Where T0_x is the piecewise T0 value, if it exists, correponding to the group the TOA belongs to. If T0_x does not exist, use the global T0 vlaue. :2: (WARNING/2) Title underline too short. Finds (barycentricTOA - T0_x). Where T0_x is the piecewise T0 value, if it exists, correponding to the group the TOA belongs to. If T0_x does not exist, use the global T0 vlaue. ---------- returns time since T0 rtype astropy.quantity.Quantity
in_piece
(par)Finds which TOAs reference which piecewise binary parameter group using the group_index_array property. :2: (WARNING/2) Title underline too short. Finds which TOAs reference which piecewise binary parameter group using the group_index_array property. ---------- par str Name of piecewise parameter e.g. 'T0X_0001' or 'A1X_0001' :5: (WARNING/2) Definition list ends without a blank line; unexpected unindent. returns boolean list (length: self._t). True where TOA references a given group, False otherwise. binary piecewise parameter label str. e.g. 'T0X' or 'A1X'. rtype list
nu
()True anomaly (Ae)
omega
()orbits
()pb
()pbdot
()pbprime
()Creates a list of piecewise orbital parameters to use in calculations. It is the same dimensions as the TOAs loaded in. Each entry is the piecewise parameter value from the group it belongs to. :2: (WARNING/2) Title underline too short. Creates a list of piecewise orbital parameters to use in calculations. It is the same dimensions as the TOAs loaded in. Each entry is the piecewise parameter value from the group it belongs to. ---------- t : Quantity. TOA, not necesserily barycentered :returns: Quantity (length: t). T0X parameter to use for each TOA in calculations. :5: (ERROR/3) Unexpected indentation. Quantity (length: t). A1X parameter to use for each TOA in calculations. :6: (WARNING/2) Block quote ends without a blank line; unexpected unindent. rtype list
prtl_der
(y, x)Find the partial derivatives in binary model pdy/pdx :param y: Name of variable to be differentiated :type y: str :param x: Name of variable the derivative respect to :type x: str
search_alias
(parname)set_param_values
([valDict])Set the parameters and assign values,
setup_internal_structures
([valDict])t0
()toa_belongs_in_group
(toas)Get the piece a TOA belongs to by finding which checking upper/lower edges of each piece. :2: (WARNING/2) Title underline too short. Get the piece a TOA belongs to by finding which checking upper/lower edges of each piece. ---------- toas : Astropy.quantity.Quantity. :returns: int (length: t). Group numbers :rtype: list
update_input
(**updates)Update the toas and parameters.
Attributes
T0
t
tt0
- 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
- piecewise_parameter_from_information_array(t)[source]
Creates a list of piecewise orbital parameters to use in calculations. It is the same dimensions as the TOAs loaded in. Each entry is the piecewise parameter value from the group it belongs to.
t : Quantity. TOA, not necesserily barycentered :returns: Quantity (length: t). T0X parameter to use for each TOA in calculations.
Quantity (length: t). A1X parameter to use for each TOA in calculations.
- rtype:
list
- toa_belongs_in_group(toas)[source]
Get the piece a TOA belongs to by finding which checking upper/lower edges of each piece.
toas : Astropy.quantity.Quantity. :returns: int (length: t). Group numbers :rtype: list
- get_group_boundaries()[source]
Get the piecewise group boundaries from the dictionary of piecewise parameter information. :returns: list (length: number of pieces). Contains all pieces’ lower edge.
list (length: number of pieces). Contains all pieces’ upper edge.
- Return type:
- get_tt0(barycentricTOA)[source]
Finds (barycentricTOA - T0_x). Where T0_x is the piecewise T0 value, if it exists, correponding to the group the TOA belongs to. If T0_x does not exist, use the global T0 vlaue.
- returns:
time since T0
- rtype:
astropy.quantity.Quantity
- in_piece(par)[source]
Finds which TOAs reference which piecewise binary parameter group using the group_index_array property.
- parstr
Name of piecewise parameter e.g. ‘T0X_0001’ or ‘A1X_0001’
- returns:
boolean list (length: self._t). True where TOA references a given group, False otherwise. binary piecewise parameter label str. e.g. ‘T0X’ or ‘A1X’.
- rtype:
list
- d_E_d_T0X()[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
- prtl_der(y, x)[source]
Find the partial derivatives in binary model pdy/pdx :param y: Name of variable to be differentiated :type y: str :param x: Name of variable the derivative respect to :type x: str
- Returns:
The derivatives pdy/pdx
- Return type:
np.array
- BTdelay()
Full BT model delay
- E()
Eccentric Anomaly
- M()
Orbit phase.
- add_binary_params(parameter, defaultValue, unit=False)
Add one parameter to the binary class.
- binary_delay()
Returns total pulsar binary delay.
- Returns:
Pulsar binary delay in the units of second
- Return type:
- compute_eccentric_anomaly(eccentricity, mean_anomaly)
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
- d_E_d_T0()
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)
derivative for E respect to binary parameter.
- Parameters:
par (string) – parameter name
- Return type:
Derivative of E respect to par
- d_M_d_par(par)[source]
derivative for M respect to bianry parameter. :param par: parameter name :type par: string
- Return type:
Derivitve of M respect to par
- d_binarydelay_d_par(par)
Get the binary delay derivatives respect to parameters.
- Parameters:
par (str) – Parameter name.
- d_nu_d_ECC()
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_T0()
Analytic calculation.
dnu/dT0 = dnu/de*de/dT0+dnu/dE*dE/dT0 de/dT0 = -EDOT
- d_nu_d_par(par)
derivative for nu respect to binary parameter.
- Parameters:
par (string) – parameter name
- Return type:
Derivative of nu respect to par
- 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.
- Parameters:
par (string) – parameter name
- Return type:
Derivative of omega respect to par
- d_pb_d_par(par)
derivative for pbprime respect to binary parameter.
- Parameters:
par (string) – parameter name
- Return type:
Derivative of M respect to par
- delayL1()
First term of Blandford & Teukolsky (1976), ApJ, 205, 580-591, eq 2.33/ First left-hand term of W.M. Smart, (1962), “Spherical Astronomy”, p35 delay equation.
alpha * (cosE-e) alpha = a1*sin(omega) Here a1 is in the unit of light second as distance
- delayL2()
Second term of Blandford & Teukolsky (1976), ApJ, 205, 580-591, eq 2.33/ / Second left-hand term of W.M. Smart, (1962), “Spherical Astronomy”, p35 delay equation.
(beta + gamma) * sinE beta = (1-e^2)*a1*cos(omega) Here a1 is in the unit of light second as distance
- delayR()
Third term of Blandford & Teukolsky (1976), ApJ, 205, 580-591, eq 2.33 / Right-hand term of W.M. Smart, (1962), “Spherical Astronomy”, p35 delay equation.
(alpha*(cosE-e)+(beta+gamma)*sinE)*(1-alpha*sinE - beta*sinE)/ (pb*(1-e*coeE)) (alpha*(cosE-e)+(beta+gamma)*sinE) is define in delayL1 and delayL2
- ecc()
Calculate eccentricity with EDOT
- nu()
True anomaly (Ae)
- update_input(**updates)
Update the toas and parameters.