enterprise.signals package

Submodules

enterprise.signals.anis_coefficients module

enterprise.signals.deterministic_signals module

enterprise.signals.gp_bases module

Utilities module containing various useful functions for use in other modules.

enterprise.signals.gp_bases.createfourierdesignmatrix_dm(toas, freqs, nmodes=30, Tspan=None, pshift=False, fref=1400, logf=False, fmin=None, fmax=None, modes=None)[source]

Construct DM-variation fourier design matrix. Current normalization expresses DM signal as a deviation [seconds] at fref [MHz]

Parameters
  • toas – vector of time series in seconds

  • freqs – radio frequencies of observations [MHz]

  • nmodes – number of fourier coefficients to use

  • Tspan – option to some other Tspan

  • pshift – option to add random phase shift

  • fref – reference frequency [MHz]

  • logf – use log frequency spacing

  • fmin – lower sampling frequency

  • fmax – upper sampling frequency

  • modes – option to provide explicit list or array of sampling frequencies

Returns

F: DM-variation fourier design matrix

Returns

f: Sampling frequencies

enterprise.signals.gp_bases.createfourierdesignmatrix_env(toas, log10_Amp=- 7, log10_Q=2.4771212547196626, t0=4579200000, nmodes=30, Tspan=None, logf=False, fmin=None, fmax=None, modes=None)[source]

Construct fourier design matrix with gaussian envelope.

Parameters
  • toas – vector of time series in seconds

  • nmodes – number of fourier coefficients to use

  • freqs – radio frequencies of observations [MHz]

  • freq – option to output frequencies

  • Tspan – option to some other Tspan

  • logf – use log frequency spacing

  • fmin – lower sampling frequency

  • fmax – upper sampling frequency

  • log10_Amp – log10 of the Amplitude [s]

  • t0 – mean of gaussian envelope [s]

  • log10_Q – log10 of standard deviation of gaussian envelope [days]

  • modes – option to provide explicit list or array of sampling frequencies

Returns

F: fourier design matrix with gaussian envelope

Returns

f: Sampling frequencies

enterprise.signals.gp_bases.createfourierdesignmatrix_eph(t, nmodes, phi, theta, freq=False, Tspan=None, logf=False, fmin=None, fmax=None, modes=None)[source]
enterprise.signals.gp_bases.createfourierdesignmatrix_ephem(toas, pos, nmodes=30, Tspan=None)[source]

Construct ephemeris perturbation Fourier design matrix and frequencies. The matrix contains nmodes*6 columns, ordered as by frequency first, Cartesian coordinate second:

sin(f0) [x], sin(f0) [y], sin(f0) [z], cos(f0) [x], cos(f0) [y], cos(f0) [z], sin(f1) [x], sin(f1) [y], sin(f1) [z], …

The corresponding frequency vector repeats every entry six times. This design matrix should be used with monopole_orf and with a powerlaw that specifies components=6.

Parameters
  • toas – vector of time series in seconds

  • pos – pulsar position as Cartesian vector

  • nmodes – number of Fourier coefficients

  • Tspan – Tspan used to define Fourier bins

Returns

F: Fourier design matrix of shape (len(toas),6*nmodes)

Returns

f: Sampling frequencies (6*nmodes)

enterprise.signals.gp_bases.createfourierdesignmatrix_red(toas, nmodes=30, Tspan=None, logf=False, fmin=None, fmax=None, pshift=False, modes=None, pseed=None)[source]

Construct fourier design matrix from eq 11 of Lentati et al, 2013 :param toas: vector of time series in seconds :param nmodes: number of fourier coefficients to use :param freq: option to output frequencies :param Tspan: option to some other Tspan :param logf: use log frequency spacing :param fmin: lower sampling frequency :param fmax: upper sampling frequency :param pshift: option to add random phase shift :param pseed: option to provide phase shift seed :param modes: option to provide explicit list or array of

sampling frequencies

Returns

F: fourier design matrix

Returns

f: Sampling frequencies

enterprise.signals.gp_priors module

Utilities module containing various useful functions for use in other modules.

enterprise.signals.gp_priors.InvGamma(alpha=1, gamma=1, size=None)[source]

Class factory for Inverse Gamma parameters.

enterprise.signals.gp_priors.InvGammaPrior(value, alpha=1, gamma=1)[source]

Prior function for InvGamma parameters.

enterprise.signals.gp_priors.InvGammaSampler(alpha=1, gamma=1, size=None)[source]

Sampling function for Uniform parameters.

enterprise.signals.gp_priors.broken_powerlaw(f, log10_A, gamma, delta, log10_fb, kappa=0.1)[source]

Generic broken powerlaw spectrum. :param f: sampling frequencies :param A: characteristic strain amplitude [set for gamma at f=1/yr] :param gamma: negative slope of PSD for f > f_break [set for comparison

at f=1/yr (default 13/3)]

Parameters
  • delta – slope for frequencies < f_break

  • log10_fb – log10 transition frequency at which slope switches from gamma to delta

  • kappa – smoothness of transition (Default = 0.1)

enterprise.signals.gp_priors.free_spectrum(f, log10_rho=None)[source]

Free spectral model. PSD amplitude at each frequency is a free parameter. Model is parameterized by S(f_i) =

ho_i^2 * T,

where

ho_i is the free parameter and T is the observation

length.

enterprise.signals.gp_priors.infinitepower(f)[source]
enterprise.signals.gp_priors.powerlaw(f, log10_A=- 16, gamma=5, components=2)[source]
enterprise.signals.gp_priors.powerlaw_genmodes(f, log10_A=- 16, gamma=5, components=2, wgts=None)[source]
enterprise.signals.gp_priors.t_process(f, log10_A=- 15, gamma=4.33, alphas=None)[source]

t-process model. PSD amplitude at each frequency is a fuzzy power-law.

enterprise.signals.gp_priors.t_process_adapt(f, log10_A=- 15, gamma=4.33, alphas_adapt=None, nfreq=None)[source]

t-process model. PSD amplitude at each frequency is a fuzzy power-law.

enterprise.signals.gp_priors.turnover(f, log10_A=- 15, gamma=4.33, lf0=- 8.5, kappa=3.3333333333333335, beta=0.5)[source]
enterprise.signals.gp_priors.turnover_knee(f, log10_A, gamma, lfb, lfk, kappa, delta)[source]

Generic turnover spectrum with a high-frequency knee. :param f: sampling frequencies of GWB :param A: characteristic strain amplitude at f=1/yr :param gamma: negative slope of PSD around f=1/yr (usually 13/3) :param lfb: log10 transition frequency at which environment dominates GWs :param lfk: log10 knee frequency due to population finiteness :param kappa: smoothness of turnover (10/3 for 3-body stellar scattering) :param delta: slope at higher frequencies

enterprise.signals.gp_signals module

Contains class factories for Gaussian Process (GP) signals. GP signals are defined as the class of signals that have a basis function matrix and basis prior vector..

enterprise.signals.gp_signals.BasisCommonGP(priorFunction, basisFunction, orfFunction, coefficients=False, combine=True, name='')[source]
enterprise.signals.gp_signals.BasisGP(priorFunction, basisFunction, coefficients=False, combine=True, selection=<class 'enterprise.signals.selections.Selection.<locals>.Selection'>, name='')[source]

Class factory for generic GPs with a basis matrix.

enterprise.signals.gp_signals.EcorrBasisModel(log10_ecorr=<class 'enterprise.signals.parameter.Uniform.<locals>.Uniform'>, coefficients=False, selection=<class 'enterprise.signals.selections.Selection.<locals>.Selection'>, name='basis_ecorr')[source]

Convenience function to return a BasisGP class with a quantized ECORR basis.

enterprise.signals.gp_signals.FourierBasisCommonGP(spectrum, orf, coefficients=False, combine=True, components=20, Tspan=None, modes=None, name='common_fourier', pshift=False, pseed=None)[source]
enterprise.signals.gp_signals.FourierBasisCommonGP_ephem(spectrum, components, Tspan, name='ephem_gp')[source]
enterprise.signals.gp_signals.FourierBasisCommonGP_physicalephem(frame_drift_rate=1e-09, d_jupiter_mass=1.5497669e-11, d_saturn_mass=8.17306184e-12, d_uranus_mass=5.71923361e-11, d_neptune_mass=7.96103855e-11, jup_orb_elements=0.05, sat_orb_elements=0.5, model='setIII', coefficients=False, name='phys_ephem_gp')[source]

Class factory for physical ephemeris corrections as a common GP. Individual perturbations can be excluded by setting the corresponding prior sigma to None.

Parameters
  • frame_drift_rate – Gaussian sigma for frame drift rate

  • d_jupiter_mass – Gaussian sigma for Jupiter mass perturbation

  • d_saturn_mass – Gaussian sigma for Saturn mass perturbation

  • d_uranus_mass – Gaussian sigma for Uranus mass perturbation

  • d_neptune_mass – Gaussian sigma for Neptune mass perturbation

  • jup_orb_elements – Gaussian sigma for Jupiter orbital elem. perturb.

  • sat_orb_elements – Gaussian sigma for Saturn orbital elem. perturb.

  • model – vector basis used by Jupiter and Saturn perturb.; see PhysicalEphemerisSignal, defaults to “setIII”

  • coefficients – if True, treat GP coefficients as enterprise parameters; if False, marginalize over them

Returns

BasisCommonGP representing ephemeris perturbations

enterprise.signals.gp_signals.FourierBasisGP(spectrum, coefficients=False, combine=True, components=20, selection=<class 'enterprise.signals.selections.Selection.<locals>.Selection'>, Tspan=None, modes=None, name='red_noise', pshift=False, pseed=None)[source]

Convenience function to return a BasisGP class with a fourier basis.

enterprise.signals.gp_signals.TimingModel(coefficients=False, name='linear_timing_model', use_svd=False, normed=True)[source]

Class factory for marginalized linear timing model signals.

enterprise.signals.gp_signals.WidebandTimingModel(dmefac=<class 'enterprise.signals.parameter.Uniform.<locals>.Uniform'>, log10_dmequad=<class 'enterprise.signals.parameter.Uniform.<locals>.Uniform'>, dmjump=<class 'enterprise.signals.parameter.Uniform.<locals>.Uniform'>, dmefac_selection=<class 'enterprise.signals.selections.Selection.<locals>.Selection'>, log10_dmequad_selection=<class 'enterprise.signals.selections.Selection.<locals>.Selection'>, dmjump_selection=<class 'enterprise.signals.selections.Selection.<locals>.Selection'>, dmjump_ref=None, name='wideband_timing_model')[source]

Class factory for marginalized linear timing model signals that take wideband TOAs and DMs. Currently assumes DMX for DM model.

enterprise.signals.gp_signals.ecorr_basis_prior(weights, log10_ecorr=- 8)[source]

Returns the ecorr prior. :param weights: A vector or weights for the ecorr prior.

enterprise.signals.parameter module

Contains parameter types for use in enterprise Signal classes.

class enterprise.signals.parameter.ConstantParameter(name)[source]

Bases: object

Constant Parameter base class.

property value
class enterprise.signals.parameter.FunctionBase[source]

Bases: object

class enterprise.signals.parameter.Parameter(name)[source]

Bases: object

get_logpdf(value=None, **kwargs)[source]
get_pdf(value=None, **kwargs)[source]
sample(**kwargs)[source]
property params
property size
enterprise.signals.parameter.Constant(val=None)[source]

Class factory for Constant parameters. Leave val=None to set value later, for example with signal_base.PTA.set_default_params().

enterprise.signals.parameter.Function(func, name='', **func_kwargs)[source]
enterprise.signals.parameter.GPCoefficients(logprior, size)[source]

Class factory for GP coefficients, which are usually created inside gp_signals.BasisGP.

enterprise.signals.parameter.LinearExp(pmin, pmax, size=None)[source]

Class factory for LinearExp parameters (with pdf(x) ~ 10^x, and 0 outside [pmin,``max``]). Handles vectors correctly if pmin and pmax are scalars or if size == len(pmin) == len(pmax)

Parameters
  • pmin – minimum of range

  • pmax – maximum of range

  • size – length for vector parameter (default None)

Returns

LinearExp parameter class

enterprise.signals.parameter.LinearExpPrior(value, pmin, pmax)[source]

Prior function for LinearExp parameters.

enterprise.signals.parameter.LinearExpSampler(pmin, pmax, size)[source]

Sampling function for LinearExp parameters.

enterprise.signals.parameter.Normal(mu=0, sigma=1, size=None)[source]

Class factory for Normal parameters (with pdf(x) ~ N(mu,``sigma``)). Handles vectors correctly if size == len(mu) == len(sigma), in which case sigma is taken as the sqrt of the diagonal of the covariance matrix; sigma can also be given passed as the size x size covariance matrix.

Parameters
  • mu – center of normal distribution

  • sigma – standard deviation of normal distribution

  • size – length for vector parameter

Returns

Normal parameter class

enterprise.signals.parameter.NormalPrior(value, mu, sigma)[source]

Prior function for Normal parameters.

enterprise.signals.parameter.NormalSampler(mu, sigma, size=None)[source]

Sampling function for Normal parameters.

enterprise.signals.parameter.Uniform(pmin, pmax, size=None)[source]

Class factory for Uniform parameters (with pdf(x) ~ 1/[pmax - pmin] inside [pmin,pmax], 0 outside. Handles vectors correctly, if pmin and pmax are scalars, or if len(size) == len(pmin) == len(pmax)

Parameters
  • pmin – minimum of uniform range

  • pmax – maximum of uniform range

  • size – length for vector parameter

Returns

Uniform parameter class

enterprise.signals.parameter.UniformPrior(value, pmin, pmax)[source]

Prior function for Uniform parameters.

enterprise.signals.parameter.UniformSampler(pmin, pmax, size=None)[source]

Sampling function for Uniform parameters.

enterprise.signals.parameter.UserParameter(prior=None, logprior=None, sampler=None, size=None)[source]

Class factory for UserParameter, implementing Enterprise parameters with arbitrary priors. The prior is specified by way of an Enterprise Function of the form prior(value, [par1, par2]). Optionally, include sampler (a function with the same parameters as prior), to allow random sampling of the parameter through enterprise.signals.parameter.sample.

Parameters
  • prior – parameter prior pdf, given as Enterprise Function

  • sampler – function returning a randomly sampled parameter according to prior

  • size – length for vector parameter

Returns

UserParameter class

enterprise.signals.parameter.function(func)[source]

Decorator for Function.

enterprise.signals.parameter.get_funcargs(func)[source]

Convenience function to get args and kwargs of any function.

enterprise.signals.parameter.sample(parlist)[source]

Sample a list of Parameters consistently (i.e., keeping track of hyperparameters).

enterprise.signals.selections module

Contains various selection functions to mask parameters by backend flags, time-intervals, etc.

enterprise.signals.selections.Selection(func)[source]

Class factory for TOA selection.

enterprise.signals.selections.by_backend(backend_flags)[source]

Selection function to split by backend flags.

enterprise.signals.selections.by_band(flags)[source]

Selection function to split by PPTA frequency band under -B flag

enterprise.signals.selections.by_frontend(flags)[source]

Selection function to split by frontend under -fe flag

enterprise.signals.selections.call_me_maybe(obj)[source]

See here for description.

enterprise.signals.selections.cut_half(toas)[source]

Selection function to split by data segment

enterprise.signals.selections.nanograv_backends(backend_flags)[source]

Selection function to split by NANOGRav backend flags only.

enterprise.signals.selections.no_selection(toas)[source]

Default selection with no splitting.

enterprise.signals.selections.selection_func(func)[source]

enterprise.signals.signal_base module

Defines the signal base classes and metaclasses. All signals will then be derived from these base classes.

class enterprise.signals.signal_base.BlockMatrix(blocks, slices, nvec=0)[source]

Bases: object

solve(other, left_array=None, logdet=False)[source]
class enterprise.signals.signal_base.CommonSignal(psr)[source]

Bases: enterprise.signals.signal_base.Signal

Base class for CommonSignal objects.

classmethod get_phicross(signal1, signal2, params)[source]
get_phiinv(params)[source]

Returns inverse of the covaraince of basis amplitudes.

class enterprise.signals.signal_base.LogLikelihood(pta)[source]

Bases: object

class enterprise.signals.signal_base.MetaCollection[source]

Bases: type

Metaclass for Signal collections. Allows addition of SignalCollection classes.

class enterprise.signals.signal_base.MetaSignal[source]

Bases: type

Metaclass for Signals. Allows addition of Signal classes.

class enterprise.signals.signal_base.PTA(init, lnlikelihood=<class 'enterprise.signals.signal_base.LogLikelihood'>)[source]

Bases: object

get_TNT(params)[source]
get_TNr(params)[source]
get_basis(params={})[source]
get_delay(params={})[source]
get_lnlikelihood(params, **kwargs)[source]
get_lnprior(params)[source]
get_logsignalprior(params)[source]
get_ndiag(params={})[source]
get_phi(params, cliques=False)[source]
get_phiinv(params, logdet=False, method='cliques')[source]
get_phiinv_byfreq_cliques(params, logdet=False, cholesky=False)[source]
get_phiinv_byfreq_partition(params, logdet=False)[source]
get_phiinv_sparse(params, logdet=False)[source]
get_rNr_logdet(params)[source]
get_residuals()[source]
get_signal(name)[source]

Returns Signal instance given the signal name.

map_params(xs)[source]
set_default_params(params)[source]
summary(include_params=True, to_stdout=False)[source]

generate summary string for PTA model

Parameters
  • include_params – [bool] list all parameters for each signal

  • to_stdout – [bool] print summary to stdout instead of returning it

Returns

[string]

property param_names
property params
property pulsarmodels
property pulsars
property signals

Return signal dictionary.

class enterprise.signals.signal_base.ShermanMorrison(jvec, slices, nvec=0.0)[source]

Bases: object

Custom container class for Sherman-morrison array inversion.

solve(other, left_array=None, logdet=False)[source]
class enterprise.signals.signal_base.Signal(psr)[source]

Bases: object

Base class for Signal objects.

get(parname, params={})[source]
get_basis(params=None)[source]

Returns the basis array of shape N_toa x N_basis.

get_delay(params)[source]

Returns the waveform of a deterministic signal.

get_logsignalprior(params)[source]

Returns an additional prior/likelihood terms associated with a signal.

get_ndiag(params)[source]

Returns the diagonal of the white noise vector N.

This method also supports block diagonal sparse matrices.

get_phi(params)[source]

Returns a diagonal covariance matrix of the basis amplitudes.

get_phiinv(params)[source]

Returns inverse of the covaraince of basis amplitudes.

set_default_params(params)[source]

Set default parameters.

property param_names
property params
class enterprise.signals.signal_base.cholesky(x)[source]

Bases: object

inv()[source]
logdet()[source]
class enterprise.signals.signal_base.csc_matrix_alt(arg1, shape=None, dtype=None, copy=False)[source]

Bases: scipy.sparse.csc.csc_matrix

Sub-class of scipy.sparse.csc_matrix with custom add and solve methods.

solve(other, left_array=None, logdet=False)[source]
class enterprise.signals.signal_base.ndarray_alt(inputarr)[source]

Bases: numpy.ndarray

Sub-class of np.ndarray with custom solve method.

solve(other, left_array=None, logdet=False)[source]
enterprise.signals.signal_base.SignalCollection(metasignals)[source]

Class factory for SignalCollection objects.

enterprise.signals.signal_base.cache_call(attrs, limit=2)[source]

This decorator caches the output of a class method that takes a single parameter ‘params’. It saves the cache in the instance attributes _cache_<methodname> and _cache_list_<methodname>.

The cache keys are listed in the class attribute (or attributes) specified in the initial decorator call. For instance, if the decorator is applied as @cache_call(‘basis_params’), then the parameters listed in self.basis_params (together with their values) will be used as the key.

The parameter ‘limit’ specifies the number of entries saved in the cache.

enterprise.signals.utils module

Utilities module containing various useful functions for use in other modules.

class enterprise.signals.utils.KernelMatrix(init)[source]

Bases: numpy.ndarray

add(other, idx)[source]
inv(logdet=False)[source]
set(other, idx)[source]
enterprise.signals.utils.anis_orf(pos1, pos2, params, **kwargs)[source]

Anisotropic GWB spatial correlation function.

enterprise.signals.utils.bwm_delay(toas, pos, log10_h=- 14.0, cos_gwtheta=0.0, gwphi=0.0, gwpol=0.0, t0=55000, antenna_pattern_fn=None)[source]

Function that calculates the earth-term gravitational-wave burst-with-memory signal, as described in: Seto et al, van haasteren and Levin, phsirkov et al, Cordes and Jenet. This version uses the F+/Fx polarization modes, as verified with the Continuous Wave and Anisotropy papers.

Parameters
  • toas – Time-of-arrival measurements [s]

  • pos – Unit vector from Earth to pulsar

  • log10_h – log10 of GW strain

  • cos_gwtheta – Cosine of GW polar angle

  • gwphi – GW azimuthal polar angle [rad]

  • gwpol – GW polarization angle

  • t0 – Burst central time [day]

  • antenna_pattern_fn – User defined function that takes pos, gwtheta, gwphi as arguments and returns (fplus, fcross)

Returns

the waveform as induced timing residuals (seconds)

enterprise.signals.utils.calculate_splus_scross(nmax, mc, dl, h0, F, e, t, l0, gamma, gammadot, inc)[source]

Calculate splus and scross for a CGW summed over all harmonics. This waveform differs slightly from that in Taylor et al (2016) in that it includes the time dependence of the advance of periastron.

Parameters
  • nmax – Total number of harmonics to use

  • mc – Chirp mass of binary [Solar Mass]

  • dl – Luminosity distance [Mpc]

  • F – Orbital frequency of binary [Hz]

  • e – Orbital Eccentricity

  • t – TOAs [s]

  • l0 – Initial eccentric anomoly [rad]

  • gamma – Angle of periastron advance [rad]

  • gammadot – Time derivative of angle of periastron advance [rad/s]

  • inc – Inclination angle [rad]

Return splus, scross

plus and cross time-domain waveforms for a CGW

enterprise.signals.utils.create_gw_antenna_pattern(pos, gwtheta, gwphi)[source]

Function to create pulsar antenna pattern functions as defined in Ellis, Siemens, and Creighton (2012). :param pos: Unit vector from Earth to pulsar :param gwtheta: GW polar angle in radians :param gwphi: GW azimuthal angle in radians

Returns

(fplus, fcross, cosMu), where fplus and fcross are the plus and cross antenna pattern functions and cosMu is the cosine of the angle between the pulsar and the GW source.

enterprise.signals.utils.create_quantization_matrix(toas, dt=1, nmin=2)[source]

Create quantization matrix mapping TOAs to observing epochs.

enterprise.signals.utils.create_stabletimingdesignmatrix(designmat, fastDesign=True)[source]

Stabilize the timing-model design matrix.

Parameters
  • designmat – Pulsar timing model design matrix

  • fastDesign – Stabilize the design matrix the fast way [True]

Returns

Mm: Stabilized timing model design matrix

enterprise.signals.utils.createfourierdesignmatrix_physicalephem(toas, planetssb, pos_t, frame_drift_rate=1e-09, d_jupiter_mass=1.5497669e-11, d_saturn_mass=8.17306184e-12, d_uranus_mass=5.71923361e-11, d_neptune_mass=7.96103855e-11, jup_orb_elements=0.05, sat_orb_elements=0.5, model='setIII')[source]

Construct physical ephemeris perturbation design matrix and ‘frequencies’. Parameters can be excluded by setting the corresponding prior sigma to None

Parameters
  • toas – vector of time series in seconds

  • pos – pulsar position as Cartesian vector

  • frame_drift_rate – normal sigma for frame drift rate

  • d_jupiter_mass – normal sigma for Jupiter mass perturbation

  • d_saturn_mass – normal sigma for Saturn mass perturbation

  • d_uranus_mass – normal sigma for Uranus mass perturbation

  • d_neptune_mass – normal sigma for Neptune mass perturbation

  • jup_orb_elements – normal sigma for Jupiter orbital elem. perturb.

  • sat_orb_elements – normal sigma for Saturn orbital elem. perturb.

  • model – vector basis used by Jupiter and Saturn perturb.; see PhysicalEphemerisSignal, defaults to “setIII”

Returns

F: Fourier design matrix of shape (len(toas), nvecs)

Returns

sigmas: Phi sigmas (nvecs, to be passed to physicalephem_spectrum)

enterprise.signals.utils.dipole_orf(pos1, pos2)[source]

Dipole spatial correlation function.

enterprise.signals.utils.dmass(planet, dm_over_Msun)[source]
enterprise.signals.utils.ecl2eq_vec(x)[source]

Rotate (n,3) vector time series from ecliptic to equatorial.

enterprise.signals.utils.eq2ecl_vec(x)[source]

Rotate (n,3) vector time series from equatorial to ecliptic.

enterprise.signals.utils.euler_vec(z, y, x, n)[source]

Return (n,3,3) tensor with each (3,3) block containing an Euler rotation with angles z, y, x. Optionally each of z, y, x can be a vector of length n.

enterprise.signals.utils.get_Fdot(F, mc, e)[source]

Compute frequency derivative from Taylor et al. (2016)

Parameters
  • F – Orbital frequency [Hz]

  • mc – Chirp mass of binary [Solar Mass]

  • e – Eccentricity of binary

Returns

dF/dt

enterprise.signals.utils.get_an(n, mc, dl, h0, F, e)[source]

Compute a_n from Eq. 22 of Taylor et al. (2016).

Parameters
  • n – Harmonic number

  • mc – Chirp mass of binary [Solar Mass]

  • dl – Luminosity distance [Mpc]

  • F – Orbital frequency of binary [Hz]

  • e – Orbital Eccentricity

Returns

a_n

enterprise.signals.utils.get_bn(n, mc, dl, h0, F, e)[source]

Compute b_n from Eq. 22 of Taylor et al. (2015).

Parameters
  • n – Harmonic number

  • mc – Chirp mass of binary [Solar Mass]

  • dl – Luminosity distance [Mpc]

  • F – Orbital frequency of binary [Hz]

  • e – Orbital Eccentricity

Returns

b_n

enterprise.signals.utils.get_cn(n, mc, dl, h0, F, e)[source]

Compute c_n from Eq. 22 of Taylor et al. (2016).

Parameters
  • n – Harmonic number

  • mc – Chirp mass of binary [Solar Mass]

  • dl – Luminosity distance [Mpc]

  • F – Orbital frequency of binary [Hz]

  • e – Orbital Eccentricity

Returns

c_n

enterprise.signals.utils.get_coefficients(pta, params, n=1, phiinv_method='cliques', common_sparse=False)[source]
enterprise.signals.utils.get_coupled_constecc_eqns(y, t, mc, e0)[source]

Computes the coupled system of differential equations from Peters (1964) and Barack & Cutler (2004). This is a system of three variables:

F: Orbital frequency [Hz] phase0: Orbital phase [rad]

Parameters
  • y – Vector of input parameters [F, e, gamma]

  • t – Time [s]

  • mc – Chirp mass of binary [Solar Mass]

Returns

array of derivatives [dF/dt, dphase/dt]

enterprise.signals.utils.get_coupled_ecc_eqns(y, t, mc, q)[source]

Computes the coupled system of differential equations from Peters (1964) and Barack & Cutler (2004). This is a system of three variables:

F: Orbital frequency [Hz] e: Orbital eccentricity gamma: Angle of precession of periastron [rad] phase0: Orbital phase [rad]

Parameters
  • y – Vector of input parameters [F, e, gamma]

  • t – Time [s]

  • mc – Chirp mass of binary [Solar Mass]

  • q – Mass ratio of binary

Returns

array of derivatives [dF/dt, de/dt, dgamma/dt, dphase/dt]

enterprise.signals.utils.get_edot(F, mc, e)[source]

Compute eccentricity derivative from Taylor et al. (2016)

Parameters
  • F – Orbital frequency [Hz]

  • mc – Chirp mass of binary [Solar Mass]

  • e – Eccentricity of binary

Returns

de/dt

enterprise.signals.utils.get_gammadot(F, mc, q, e)[source]

Compute gamma dot from Barack and Cutler (2004)

Parameters
  • F – Orbital frequency [Hz]

  • mc – Chirp mass of binary [Solar Mass]

  • q – Mass ratio of binary

  • e – Eccentricity of binary

Returns

dgamma/dt

enterprise.signals.utils.get_planet_orbital_elements(model='setIII')[source]

Grab physical ephemeris model files

enterprise.signals.utils.hd_orf(pos1, pos2)[source]

Hellings & Downs spatial correlation function.

enterprise.signals.utils.linear_interp_basis(toas, dt=2592000)[source]

Provides a basis for linear interpolation.

Parameters
  • toas – Pulsar TOAs in seconds

  • dt – Linear interpolation step size in seconds.

Returns

Linear interpolation basis and nodes

enterprise.signals.utils.make_ecc_interpolant()[source]

Make interpolation function from eccentricity file to determine number of harmonics to use for a given eccentricity.

Returns

interpolant

enterprise.signals.utils.monopole_orf(pos1, pos2)[source]

Monopole spatial correlation function.

enterprise.signals.utils.normed_tm_basis(Mmat, norm=None)[source]
enterprise.signals.utils.physical_ephem_delay(toas, planetssb, pos_t, frame_drift_rate=0, d_jupiter_mass=0, d_saturn_mass=0, d_uranus_mass=0, d_neptune_mass=0, jup_orb_elements=array([0.0, 0.0, 0.0, 0.0, 0.0, 0.0]), sat_orb_elements=array([0.0, 0.0, 0.0, 0.0, 0.0, 0.0]), times=None, jup_orbit=None, sat_orbit=None, equatorial=True)[source]
enterprise.signals.utils.physicalephem_spectrum(sigmas)[source]
enterprise.signals.utils.quant2ind(U)[source]

Use quantization matrix to return slices of non-zero elements.

Parameters

U – quantization matrix

Returns

list of `slice`s for non-zero elements of U

Note

This function assumes that the pulsar TOAs were sorted by time.

enterprise.signals.utils.solve_coupled_constecc_solution(F0, e0, phase0, mc, t)[source]

Compute the solution to the coupled system of equations from from Peters (1964) and Barack & Cutler (2004) at a given time.

Parameters
  • F0 – Initial orbital frequency [Hz]

  • mc – Chirp mass of binary [Solar Mass]

  • t – Time at which to evaluate solution [s]

Returns

(F(t), phase(t))

enterprise.signals.utils.solve_coupled_ecc_solution(F0, e0, gamma0, phase0, mc, q, t)[source]

Compute the solution to the coupled system of equations from from Peters (1964) and Barack & Cutler (2004) at a given time.

Parameters
  • F0 – Initial orbital frequency [Hz]

  • e0 – Initial orbital eccentricity

  • gamma0 – Initial angle of precession of periastron [rad]

  • mc – Chirp mass of binary [Solar Mass]

  • q – Mass ratio of binary

  • t – Time at which to evaluate solution [s]

Returns

(F(t), e(t), gamma(t), phase(t))

enterprise.signals.utils.ss_framerotate(mjd, planet, x, y, z, dz, offset=None, equatorial=False)[source]

Rotate planet trajectory given as (n,3) tensor, by ecliptic Euler angles x, y, z, and by z rate dz. The rate has units of rad/year, and is referred to offset 2010/1/1. dates must be given in MJD.

enterprise.signals.utils.svd_tm_basis(Mmat)[source]
enterprise.signals.utils.tm_prior(weights)[source]
enterprise.signals.utils.unnormed_tm_basis(Mmat)[source]

enterprise.signals.white_signals module

Contains class factories for white noise signals. White noise signals are defined as the class of signals that only modifies the white noise matrix N.

enterprise.signals.white_signals.EcorrKernelNoise(log10_ecorr=<class 'enterprise.signals.parameter.Uniform.<locals>.Uniform'>, selection=<class 'enterprise.signals.selections.Selection.<locals>.Selection'>, method='sherman-morrison', name='')[source]

Class factory for ECORR type noise.

Parameters
  • log10_ecorrParameter type for log10 or ecorr parameter.

  • selectionSelection object specifying masks for backends, time segments, etc.

  • method – Method for computing noise covariance matrix. Options include sherman-morrison, sparse, and block

Returns

EcorrKernelNoise class.

ECORR is a noise signal that is used for data with multi-channel TOAs that are nearly simultaneous in time. It is a white noise signal that is uncorrelated epoch to epoch but completely correlated for TOAs in a given observing epoch.

For this implementation we use this covariance matrix as part of the white noise covariance matrix \(N\). It can be seen from above that this covariance is block diagonal, thus allowing us to exploit special methods to make matrix manipulations easier.

In this signal implementation we offer three methods of performing these matrix operations:

sherman-morrison

Uses the Sherman-Morrison forumla to compute the matrix inverse and other matrix operations. Note: This method can only be used for covariances that can be constructed by the outer product of two vectors, \(uv^T\).

sparse

Uses Scipy Sparse matrices to construct the block diagonal covariance matrix and perform matrix operations.

block

Uses a custom scheme that uses the individual blocks from the block diagonal matrix to perform fast matrix inverse and other solve operations.

Note

The sherman-morrison method is the fastest, followed by the block and then sparse methods, however; the block and sparse methods are more general and should be used if sub-classing this signal for more complicated blocks.

enterprise.signals.white_signals.EquadNoise(log10_equad=<class 'enterprise.signals.parameter.Uniform.<locals>.Uniform'>, selection=<class 'enterprise.signals.selections.Selection.<locals>.Selection'>, name='')[source]

Class factory for EQUAD type measurement noise.

enterprise.signals.white_signals.MeasurementNoise(efac=<class 'enterprise.signals.parameter.Uniform.<locals>.Uniform'>, selection=<class 'enterprise.signals.selections.Selection.<locals>.Selection'>, name='')[source]

Class factory for EFAC type measurement noise.

enterprise.signals.white_signals.WhiteNoise(varianceFunction, selection=<class 'enterprise.signals.selections.Selection.<locals>.Selection'>, name='')[source]

Class factory for generic white noise signals.

enterprise.signals.white_signals.efac_ndiag(toaerrs, efac=1.0)[source]
enterprise.signals.white_signals.equad_ndiag(toas, log10_equad=- 8)[source]