pint.models.solar_wind_dispersion.SolarWindDispersionX

class pint.models.solar_wind_dispersion.SolarWindDispersionX[source]

Bases: SolarWindDispersionBase

This class provides a SWX model - multiple Solar Wind segments.

This model lets the user specify time ranges and fit for a different SWXDM (max solar wind DM) value and SWXP (radial power-law index) in each time range.

The default radial power-law index value of 2 corresponds to the Edwards et al. model. Other values are for the You et al./Hazboun et al. model.

Note that unlike the standard model, this model goes to 0 at opposition (the minimum). So it really represents a Delta DM. This is to make it easier to join multiple segments. However, to get the peak to still be the requested max DM the values are scaled compared to the standard model: the standard model goes from opposition (min) to conjunction (max), while this model goes from 0 to conjunction (max), so the scaling is ((conjunction - opposition)/conjuction).

See Solar Wind Examples.

To compare against a standard model:

Example

>>> swxmodel.SWXDM_0001.quantity = model.get_max_dm() * swxmodel.get_swscalings()[0]
>>> np.allclose(swxmodel.get_min_dms(),model.get_min_dm()*swxmodel.get_swscalings())
>>> toas = pint.simulation.make_fake_toas_uniform(54000, 54000 + 2 * 365, 500, model=model, obs="gbt")
>>> np.allclose((swxmodel.swx_dm(toas) + model.get_min_dm()), model.solar_wind_dm(toas))

Parameters supported:

Name / Aliases

Description

Kind

SWXDM_{number}

Max Solar Wind DM

pc / cm3

SWXP_{number}

Solar wind power-law index

number

SWXR1_{number}

Beginning of SWX interval

d

SWXR2_{number}

End of SWX interval

d

References

  • Edwards et al. 2006, MNRAS, 372, 1549; Setion 2.5.4

  • Madison et al. 2019, ApJ, 872, 150; Section 3.1.

  • Hazboun et al. (2022, ApJ, 929, 39)

  • You et al. (2012, MNRAS, 422, 1160)

Methods

add_param(param[, deriv_func, setup])

Add a parameter to the Component.

add_swx_range(mjd_start, mjd_end[, index, ...])

Add SWX range to a dispersion model with specified start/end MJD, SWXDM, and power-law index

conjunction_solar_wind_geometry(p)

Return the geometry (integrated path-length) of solar wind dispersion at conjunction

d_conjunction_solar_wind_geometry_d_swxp(p)

Derivative of conjunction_solar_wind_geometry (path length) wrt power-law index p

d_delay_d_dmparam(toas, param_name[, acc_delay])

Derivative of delay wrt to DM parameter.

d_delay_d_swxdm(toas, param_name[, acc_delay])

d_delay_d_swxp(toas, param_name[, acc_delay])

d_dm_d_swxdm(toas, param_name[, acc_delay])

d_dm_d_swxp(toas, param_name[, acc_delay])

d_opposition_solar_wind_geometry_d_swxp(p)

Derivative of opposition_solar_wind_geometry (path length) wrt power-law index p

d_solar_wind_geometry_d_swxp(toas, p)

Derivative of solar_wind_geometry (path length) wrt power-law index p

dispersion_slope_value(toas)

dispersion_time_delay(DM, freq)

Return the dispersion time delay for a set of frequency.

dispersion_type_delay(toas)

dm_value(toas)

Compute modeled DM value at given TOAs.

get_indices()

Returns an array of integers corresponding to SWX parameters.

get_max_dms()

Return approximate maximum DMs for each segment from the Solar Wind (at conjunction)

get_min_dms()

Return approximate minimum DMs for each segment from the Solar Wind (at opposition).

get_ne_sws()

Return Solar Wind electron densities at 1 AU for each segment

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.

get_swscalings()

Return the approximate scaling between the SWX model and the standard model

get_theta0()

Get elongation at conjunction

is_in_parfile(para_dict)

Check if this subclass included in parfile.

match_param_aliases(alias)

Return the parameter corresponding to this alias.

opposition_solar_wind_geometry(p)

Return the geometry (integrated path-length) of solar wind dispersion at opposition

param_help()

Print help lines for all available parameters in 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.

register_dm_deriv_funcs(func, param)

Register the derivative function in to the deriv_func dictionaries.

remove_param(param)

Remove a parameter from the Component.

remove_swx_range(index)

Removes all SWX parameters associated with a given index/list of indices.

set_ne_sws(ne_sws)

Set the DMMAXs based on an input NE_SW values (electron density at 1 AU)

set_special_params(spcl_params)

setup()

Finalize construction loaded values.

solar_wind_geometry(toas, p)

Return the geometry of solar wind dispersion (integrated path-length)

swx_delay(toas[, acc_delay])

This is a wrapper function for interacting with the TimingModel class

swx_dm(toas)

Return solar wind Delta DM for given TOAs

validate()

Validate the SWX parameters.

validate_toas(toas)

Check that this model component has TOAs where needed.

Attributes

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

theta0

solar_wind_geometry(toas, p)[source]

Return the geometry of solar wind dispersion (integrated path-length)

Implements Eqn. 11 of Hazboun et al. (2022)

Parameters:
Return type:

astropy.quantity.Quantity

conjunction_solar_wind_geometry(p)[source]

Return the geometry (integrated path-length) of solar wind dispersion at conjunction

Implements Eqn. 11 of Hazboun et al. (2022)

Parameters:

p (float) – Radial power-law index

Return type:

astropy.quantity.Quantity

opposition_solar_wind_geometry(p)[source]

Return the geometry (integrated path-length) of solar wind dispersion at opposition

Implements Eqn. 11 of Hazboun et al. (2022)

Parameters:

p (float) – Radial power-law index

Return type:

astropy.quantity.Quantity

d_solar_wind_geometry_d_swxp(toas, p)[source]

Derivative of solar_wind_geometry (path length) wrt power-law index p

The evaluation is done using Eqn. 12 in Hazboun et al. (2022). The first term involving hypergeometric functions (_hypergeom_function()) has the derivative computed approximately using a Pade expansion (_d_hypergeom_function_dp()). The second uses gamma functions (_gamma_function()) and has the derivative computed using polygamma functions (_d_gamma_function_dp()).

Parameters:
Return type:

astropy.quantity.Quantity

d_conjunction_solar_wind_geometry_d_swxp(p)[source]

Derivative of conjunction_solar_wind_geometry (path length) wrt power-law index p

The evaluation is done using Eqn. 12 in Hazboun et al. (2022). The first term involving hypergeometric functions (_hypergeom_function()) has the derivative computed approximately using a Pade expansion (_d_hypergeom_function_dp()). The second uses gamma functions (_gamma_function()) and has the derivative computed using polygamma functions (_d_gamma_function_dp()).

d_opposition_solar_wind_geometry_d_swxp(p)[source]

Derivative of opposition_solar_wind_geometry (path length) wrt power-law index p

The evaluation is done using Eqn. 12 in Hazboun et al. (2022). The first term involving hypergeometric functions (_hypergeom_function()) has the derivative computed approximately using a Pade expansion (_d_hypergeom_function_dp()). The second uses gamma functions (_gamma_function()) and has the derivative computed using polygamma functions (_d_gamma_function_dp()).

add_swx_range(mjd_start, mjd_end, index=None, swxdm=0, swxp=2, frozen=True)[source]

Add SWX range to a dispersion model with specified start/end MJD, SWXDM, and power-law index

Parameters:
  • mjd_start (float or astropy.quantity.Quantity or astropy.time.Time) – MJD for beginning of DMX event.

  • mjd_end (float or astropy.quantity.Quantity or astropy.time.Time) – MJD for end of DMX event.

  • index (int, None) – Integer label for DMX event. If None, will increment largest used index by 1.

  • swxdm (float or astropy.quantity.Quantity) – Max solar wind DM

  • swxp (float or astropy.quantity.Quantity) – Solar wind power-law index

  • frozen (bool) – Indicates whether SWXDM and SWXP will be fit.

Returns:

index – Index that has been assigned to new SWX event.

Return type:

int

remove_swx_range(index)[source]

Removes all SWX parameters associated with a given index/list of indices.

Parameters:

index (float, int, list, np.ndarray) – Number or list/array of numbers corresponding to SWX indices to be removed from model.

get_indices()[source]

Returns an array of integers corresponding to SWX parameters.

Returns:

inds – Array of SWX indices in model.

Return type:

np.ndarray

setup()[source]

Finalize construction loaded values.

validate()[source]

Validate the SWX parameters.

get_theta0()[source]

Get elongation at conjunction

Return type:

astropy.quantity.Quantity

validate_toas(toas)[source]

Check that this model component has TOAs where needed.

swx_dm(toas)[source]

Return solar wind Delta DM for given TOAs

swx_delay(toas, acc_delay=None)[source]

This is a wrapper function for interacting with the TimingModel class

print_par(format='pint')[source]
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

get_swscalings()[source]

Return the approximate scaling between the SWX model and the standard model

Return type:

np.ndarray

get_max_dms()[source]

Return approximate maximum DMs for each segment from the Solar Wind (at conjunction)

Simplified model that assumes a circular orbit

Return type:

astropy.quantity.Quantity

get_min_dms()[source]

Return approximate minimum DMs for each segment from the Solar Wind (at opposition).

Simplified model that assumes a circular orbit

Note that this value has been subtracted off of the model

Return type:

astropy.quantity.Quantity

get_ne_sws()[source]

Return Solar Wind electron densities at 1 AU for each segment

Simplified model that assumes a circular orbit

Return type:

astropy.quantity.Quantity

set_ne_sws(ne_sws)[source]

Set the DMMAXs based on an input NE_SW values (electron density at 1 AU)

Parameters:

ne_sws (astropy.quantity.Quantity) – Desired NE_SWs (should be scalar or the same length as the number of segments)

:raises ValueError : If length of input does ont match number of segments:

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.

d_delay_d_dmparam(toas, param_name, acc_delay=None)

Derivative of delay wrt to DM parameter.

Parameters:
  • toas (pint.TOAs object.) – Input toas.

  • param_name (str) – Derivative parameter name

  • acc_delay (astropy.quantity or numpy.ndarray) – Accumulated delay values. This parameter is to keep the unified API, but not used in this function.

dispersion_time_delay(DM, freq)

Return the dispersion time delay for a set of frequency.

This equation if cited from Duncan Lorimer, Michael Kramer, Handbook of Pulsar Astronomy, Second edition, Page 86, Equation [4.7] Here we assume the reference frequency is at infinity and the EM wave frequency is much larger than plasma frequency.

dm_value(toas)

Compute modeled DM value at given TOAs.

Parameters:

toas (TOAs object or TOA table(TOAs.table)) –

If given a TOAs object, it will use the whole TOA table in the

TOAs object.

Return type:

DM values at given TOAs in the unit of DM.

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.

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

register_dm_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.