This Jupyter notebook can be downloaded from noise-fitting-example.ipynb, or viewed as a python script at noise-fitting-example.py.
PINT Noise Fitting Examples
[1]:
from pint.models import get_model
from pint.simulation import make_fake_toas_uniform
from pint.logging import setup as setup_log
from pint.fitter import Fitter
import numpy as np
from io import StringIO
from astropy import units as u
from matplotlib import pyplot as plt
[2]:
setup_log(level="WARNING")
[2]:
1
Fitting for EFAC and EQUAD
[3]:
# Let us begin by simulating a dataset with an EFAC and an EQUAD.
# Note that the EFAC and the EQUAD are set as fit parameters ("1").
par = """
PSR TEST1
RAJ 05:00:00 1
DECJ 15:00:00 1
PEPOCH 55000
F0 100 1
F1 -1e-15 1
EFAC tel gbt 1.3 1
EQUAD tel gbt 1.1 1
TZRMJD 55000
TZRFRQ 1400
TZRSITE gbt
EPHEM DE440
CLOCK TT(BIPM2019)
UNITS TDB
"""
m = get_model(StringIO(par))
ntoas = 200
# EFAC and EQUAD cannot be measured separately if all TOA uncertainties
# are the same. So we must set a different toa uncertainty for each TOA.
# This is how it is in real datasets anyway.
toaerrs = np.random.uniform(0.5, 2, ntoas) * u.us
t = make_fake_toas_uniform(
startMJD=54000,
endMJD=56000,
ntoas=ntoas,
model=m,
obs="gbt",
error=toaerrs,
add_noise=True,
include_bipm=True,
)
[4]:
# Now create the fitter. The `Fitter.auto()` function creates a
# Downhill fitter. Noise parameter fitting is only available in
# Downhill fitters.
ftr = Fitter.auto(t, m)
[5]:
# Now do the fitting.
ftr.fit_toas()
[6]:
# Print the post-fit model. We can see that the EFAC and EQUAD have been
# and the uncertainties are listed.
print(ftr.model)
# Created: 2025-08-22T08:58:24.121461
# PINT_version: 1.1.4+38.g086c46d
# User: docs
# Host: build-29289377-project-85767-nanograv-pint
# OS: Linux-6.8.0-1029-aws-x86_64-with-glibc2.35
# Python: 3.11.12 (main, May 6 2025, 10:45:53) [GCC 11.4.0]
# Format: pint
# read_time: 2025-08-22T08:58:22.113329
# allow_tcb: False
# convert_tcb: False
# allow_T2: False
PSR TEST1
EPHEM DE440
CLOCK TT(BIPM2019)
UNITS TDB
START 53999.9999999862516320
FINISH 56000.0000000055922917
DILATEFREQ N
DMDATA N
NTOA 200
CHI2 199.99999962360187
CHI2R 1.036269428101564
TRES 1.8553751448691252939
RAJ 4:59:59.99999698 1 0.00000670781500043025
DECJ 14:59:59.99893229 1 0.00058393656059112968
PMRA 0.0
PMDEC 0.0
PX 0.0
F0 99.99999999999972651 1 2.6349013498996258763e-13
F1 -9.999951007540373814e-16 1 1.2045139502422000082e-20
PEPOCH 55000.0000000000000000
PLANET_SHAPIRO N
EFAC tel gbt 1.1733556269884733 1 0.18200939148321554
EQUAD tel gbt 1.11974928372613 1 0.3278883241643408
TZRMJD 55000.0000000000000000
TZRSITE gbt
TZRFRQ 1400.0
[7]:
# Let us plot the injected and measured noise parameters together to
# compare them.
plt.scatter(m.EFAC1.value, m.EQUAD1.value, label="Injected", marker="o", color="blue")
plt.errorbar(
ftr.model.EFAC1.value,
ftr.model.EQUAD1.value,
xerr=ftr.model.EFAC1.uncertainty_value,
yerr=ftr.model.EQUAD1.uncertainty_value,
marker="+",
label="Measured",
color="red",
)
plt.xlabel("EFAC_tel_gbt")
plt.ylabel("EQUAD_tel_gbt (us)")
plt.legend()
plt.show()
Fitting for ECORRs
[8]:
# Note the explicit offset (PHOFF) in the par file below.
# Implicit offset subtraction is typically not accurate enough when
# ECORR (or any other type of correlated noise) is present.
# i.e., PHOFF should be a free parameter when ECORRs are being fit.
par = """
PSR TEST2
RAJ 05:00:00 1
DECJ 15:00:00 1
PEPOCH 55000
F0 100 1
F1 -1e-15 1
PHOFF 0 1
EFAC tel gbt 1.3 1
ECORR tel gbt 1.1 1
TZRMJD 55000
TZRFRQ 1400
TZRSITE gbt
EPHEM DE440
CLOCK TT(BIPM2019)
UNITS TDB
"""
m = get_model(StringIO(par))
# ECORRs only apply when there are multiple TOAs per epoch.
# This can be simulated by providing multiple frequencies and
# setting the `multi_freqs_in_epoch` option. The `add_correlated_noise`
# option should also be set because correlated noise components
# are not simulated by default.
ntoas = 500
toaerrs = np.random.uniform(0.5, 2, ntoas) * u.us
freqs = np.linspace(1300, 1500, 4) * u.MHz
t = make_fake_toas_uniform(
startMJD=54000,
endMJD=56000,
ntoas=ntoas,
model=m,
obs="gbt",
error=toaerrs,
freq=freqs,
add_noise=True,
add_correlated_noise=True,
include_bipm=True,
multi_freqs_in_epoch=True,
)
[9]:
ftr = Fitter.auto(t, m)
[10]:
ftr.fit_toas()
[10]:
True
[11]:
print(ftr.model)
# Created: 2025-08-22T08:58:33.694027
# PINT_version: 1.1.4+38.g086c46d
# User: docs
# Host: build-29289377-project-85767-nanograv-pint
# OS: Linux-6.8.0-1029-aws-x86_64-with-glibc2.35
# Python: 3.11.12 (main, May 6 2025, 10:45:53) [GCC 11.4.0]
# Format: pint
# read_time: 2025-08-22T08:58:24.241292
# allow_tcb: False
# convert_tcb: False
# allow_T2: False
PSR TEST2
EPHEM DE440
CLOCK TT(BIPM2019)
UNITS TDB
START 53999.9999999862241435
FINISH 55984.0000000565068750
DILATEFREQ N
DMDATA N
NTOA 500
CHI2 500.0012826810599
CHI2R 1.016262769676951
TRES 1.7264700490040013856
RAJ 4:59:59.99999747 1 0.00000640581459827326
DECJ 15:00:00.00068109 1 0.00055253021032315784
PMRA 0.0
PMDEC 0.0
PX 0.0
F0 99.99999999999991027 1 2.4960591975401273373e-13
F1 -1.0000149916149782179e-15 1 1.1295046159769979103e-20
PEPOCH 55000.0000000000000000
PLANET_SHAPIRO N
EFAC tel gbt 1.2653502794719307 1 0.04632014731514219
ECORR tel gbt 1.2281354498611063 1 0.10212807063044686
TZRMJD 55000.0000000000000000
TZRSITE gbt
TZRFRQ 1400.0
PHOFF 1.364044463830738e-05 1 1.8774376602473853e-05
[12]:
# Let us plot the injected and measured noise parameters together to
# compare them.
plt.scatter(m.EFAC1.value, m.ECORR1.value, label="Injected", marker="o", color="blue")
plt.errorbar(
ftr.model.EFAC1.value,
ftr.model.ECORR1.value,
xerr=ftr.model.EFAC1.uncertainty_value,
yerr=ftr.model.ECORR1.uncertainty_value,
marker="+",
label="Measured",
color="red",
)
plt.xlabel("EFAC_tel_gbt")
plt.ylabel("ECORR_tel_gbt (us)")
plt.legend()
plt.show()