Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Disable TZRMJD corrections if TZRSITE is bat #1707

Closed
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
32 changes: 22 additions & 10 deletions src/pint/models/absolute_phase.py
Original file line number Diff line number Diff line change
Expand Up @@ -108,16 +108,28 @@ def get_TZR_toa(self, toas):
# ephem=toas.ephem,
# planets=toas.planets,
# )
tz = toa.get_TOAs_array(
(self.TZRMJD.quantity.jd1 - 2400000.5, self.TZRMJD.quantity.jd2),
obs=self.TZRSITE.value,
freqs=self.TZRFRQ.quantity,
include_bipm=clkc_info["include_bipm"],
include_gps=clkc_info["include_gps"],
ephem=toas.ephem,
planets=toas.planets,
tzr=True,
)
if self.TZRSITE.value != "@" and self.TZRSITE.value.lower() != "bat":
tz = toa.get_TOAs_array(
Comment on lines +111 to +112
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

The barycenter has more aliases than "@" and bat. It may be better to change this comparison to
self.TZRSITE.value.lower() not in [get_observatory("@" ).name] + get_observatory("@" ).aliases

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Also, please add a similar condition for the geocenter ("geo") as well.

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Hi. In tempo2's https://bitbucket.org/psrsoft/tempo2/src/master/refphs.C, this is what it does:

 21             if (obs->telID[0]=='@'            /* Have barycentric arrival time */
 22                     || strcasecmp(obs->telID,"bat")==0)
 23             {
 24                 obs->clockCorr=0;  /* therefore don't do clock corrections */
 25                 obs->delayCorr=0;
 26             }
 27             else if (strcmp(obs->telID,"STL")==0)
 28             {
 29                 obs->clockCorr=0;  /* don't do clock corrections */
 30                 obs->delayCorr=1;
 31             }
 32             else if (strcmp(obs->telID,"STL_FBAT")==0)
 33             {
 34                 obs->clockCorr=0;  /* don't do clock corrections */
 35                 obs->delayCorr=1;
 36             }
 37             else
 38             {
 39                 obs->clockCorr=1;
 40                 obs->delayCorr=1;
 41             }

It doesn't involve geocentric TOAs at all, and indeed, a par file with "COE" as TZRSITE will get the clock corrections in both PINT and tempo2. If we now disable them in PINT, I am afraid we will end up in the mirrored situation in which PINT doesn't apply the clock corrections to geocentric TZRMJD, whereas tempo2 does, producing again absolute phase differences between the two packages.

Regarding the aliases for "@" and "BAT", this can be easily changed. I have little experience with pytest, but I can try and build a test by looking at other ones in the tests folder.

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

OK... Then you can ignore the part about the geocenter.

About STL and STL_FBAT, I don't know what they are. They are not supported by PINT in any case.

Copy link
Contributor

@abhisrkckl abhisrkckl Jan 3, 2024

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I dug into this a bit more. The include_bipm and include_gps attributes of the BarycenterObs class are explicitly set to False in https://github.com/nanograv/PINT/blob/master/src/pint/observatory/special_locations.py

So this bug should not have occurred in the first place.

I think the bug is in the get_observatory() function. When I do print(get_observatory("@").include_bipm), I get True. This is because the get_observatory() function wrongly sets this from default arguments. I think the correction should be in get_observatory(). This special case treatment is not actually needed.

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Thanks. If I set defaults to None and reset them only if they are provided, as you suggest, it doesn't actually work, because in pint.toa get_observatory() gets called with explicit arguments:

2206         for obs, grp in self.get_obs_groups():
2207             site = get_observatory(
2208                 obs,
2209                 include_gps=include_gps,
2210                 include_bipm=include_bipm,
2211                 bipm_version=bipm_version,
2212             )

Any idea?

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I meant changing the get_observatory function itself to something like this:

def get_observatory(
    name, include_gps=True, include_bipm=True, bipm_version=bipm_default
):
    """Convenience function to get observatory object with options.

    This function will simply call the ``Observatory.get`` method but
    will manually modify the global observatory object after the method is called.
    Name-matching is case-insensitive.

    If the observatory is not present in the PINT list, will fallback to astropy.

    Parameters
    ----------
    name : str
        The name of the observatory
    include_gps : bool, optional
        Set False to disable UTC(GPS)->UTC clock correction.
    include_bipm : bool, optional
        Set False to disable TAI TT(BIPM) clock correction.
    bipm_version : str, optional
        Set the version of TT BIPM clock correction files.

    .. note:: This function can and should be expanded if more clock
        file switches/options are added at a public API level.

    """
    site = Observatory.get(name)
    if include_gps is not None:
        site.include_gps = include_gps
    if include_bipm is not None:
        site.include_bipm = include_bipm
        site.bipm_version = bipm_version
    return site

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Hi,
Yes. That is precisely what I did. It still fails, even if I set the defaults of get_observatory to None. This is because toa.py:get_TOAs() sets include_bipm and include_gps to True even when it is called without them (and should therefore default to None), if all TOAs provided do not have yet the clock corrections (as is the case):

 283     if all("clkcorr" not in f for f in t.table["flags"]):
 284         if include_gps is None:
 285             include_gps = True
 286         if bipm_version is None:
 287             bipm_version = bipm_default
 288         if include_bipm is None:
 289             include_bipm = True
 290         # FIXME: should we permit existing clkcorr flags?
 291         t.apply_clock_corrections(
 292             include_gps=include_gps,
 293             include_bipm=include_bipm,
 294             bipm_version=bipm_version,
 295             limits=limits,
 296         )

The conditions in lines 284 and 288 are met by default, i.e. when calling get_TOAs without explicitly setting include_gps and include_bipm, and those parameters are passed over to apply_clock_corrections, which in turn passes them to get_observatory. Therefore, even if a TOA in a tim file is already defined in the barycenter, whenever it is read with get_TOAs it will get the clock corrections, unless the user manually sets get_TOAs(..., include_gps=False, include_bipm=False). This would, however, affect all TOAs in the tim file, not only the barycentered ones!

The behavior with TZRMJD when TZRSITE is BAT follows from the same origin. While building the TZR TOA, get_TZR_toa(self,toas) in absolute_phase.py reads the clock correction information from the provided toas, and calls get_TOAs_array() explicitly, passing the include_gps and include_bipm parameters that were applied on the toas. Since these were set to True when reading them for the first time, it ends up applying the correction also on TZRMJD, even if it is defined in BAT..

At least, this was my attempt to follow the code around. It involves an interplay of get_observatory() [observatory/init.py], get_TOAs(), get_TOAs_array(), apply_clock_corrections.py [toa.py], and get_TZR_toa.py [absolute_phase.py].

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Hi @abhisrkckl ,

Unfortunately I am on a technical shift to the telescopes in La Palma, I don't have much time these days to continue digging on this. Did you have a chance to have a look at the problem?

So far, the only solution that worked consistently is adding a check for the BarycenterObs class in get_observatory() (#1707 (comment)), and it passes the tests I've implemented. but it also has some drawbacks. Namely, if apply_clock_corrections() is called on some TOAs multiple times, and there is a barycentered TOA in the list, all calls after the first one will fail.

Cheers,
Giovanni

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

No problem. I'll close this PR and open a separate one.

(self.TZRMJD.quantity.jd1 - 2400000.5, self.TZRMJD.quantity.jd2),
obs=self.TZRSITE.value,
freqs=self.TZRFRQ.quantity,
include_bipm=clkc_info["include_bipm"],
include_gps=clkc_info["include_gps"],
ephem=toas.ephem,
planets=toas.planets,
tzr=True,
)
else:
tz = toa.get_TOAs_array(
(self.TZRMJD.quantity.jd1 - 2400000.5, self.TZRMJD.quantity.jd2),
obs=self.TZRSITE.value,
freqs=self.TZRFRQ.quantity,
include_bipm=False,
include_gps=False,
ephem=toas.ephem,
planets=toas.planets,
tzr=True,
)
log.debug("Done with TZR_toa")
self.tz_cache = tz
Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@abhisrkckl I think in here self.tz_clck_info should also reflect the fact that if TZRMJD is barycentered, no corrections have been applied. What do you think?

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Sorry, it would be line 124/136, it is just outside the commit and cannot comment directly on it.

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Please see my comment above. I think if we fix the bug in get_observatory, this should get fixed automatically.

self.tz_hash = hash((self.TZRMJD.value, self.TZRSITE.value, self.TZRFRQ.value))
Expand Down
126 changes: 126 additions & 0 deletions tests/test_bary_corrections.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,126 @@
from pint.models import get_model
from pint.toa import get_TOAs
from io import StringIO
import pytest

import astropy.units as u
import numpy as np

par = """
CLK TT(BIPM2021)
EPHEM DE421
RAJ 0:00:00
DECJ 0:00:00
UNITS TDB
F0 100
PEPOCH 60000
TZRMJD 60000
TZRSITE @
TZRFRQ inf
"""

tim = """
FORMAT 1
ref 0.0 60000.0 0.0 bat
coe 0.0 60000.1 0.0 coe
foo 0.0 60000.2 0.0 pks
bar 0.0 60000.3 0.0 lst
baz 0.0 60000.3 0.0 jb
"""

bary_tempo2 = [
"60000.000000000000000",
"60000.095531198261067",
"60000.195527189099078",
"60000.295522704703600",
"60000.295522772151454",
]

phase_tempo2 = [
"0.0",
"0.552975622109443065710",
"0.913816033558759954758",
"0.168639107672788668424",
"0.751388562124702730216",
]

bary_tempo2 = u.d * np.loadtxt(bary_tempo2, dtype=np.longdouble)
phase_tempo2 = (
u.dimensionless_unscaled
* np.loadtxt(phase_tempo2, dtype=np.longdouble)
% (np.longdouble("1"))
)


@pytest.fixture
def model():
return get_model(StringIO(par))


@pytest.fixture
def toas():
return get_TOAs(StringIO(tim))


def test_first_toa_is_bat(toas):
"""Test if the ordering was respected, a condition for following tests."""
assert (
toas.__getitem__(index=0).to_TOA_list()[0].obs == "barycenter"
), "First TOA in the list should be barycentered."


def test_bat_toa_corrections(toas):
"""Test if a TOA defined in the barycenter gets the clock corrections."""
failstr = "Clock corrections applied to TOA even if in: {0}.".format(
toas.__getitem__(index=0).to_TOA_list()[0].obs
)
assert "clkcorr" not in toas.__getitem__(index=0).table["flags"][0], failstr


def test_bat_tzr_corrections(toas, model):
"""Test if a TZR defined in the barycenter gets the clock corrections."""

# Make sure the component is untouched.
assert model.components["AbsPhase"].tz_cache is None, "tz_cache already touched."

tz = model.components["AbsPhase"].get_TZR_toa(toas)
# At this point its tzr object should be True.
assert tz.tzr, f"The TZR has tzr attribute: {tzr}."

# Test wether clock corrections have been applied.
failstr = "Clock corrections applied to TZR even if in: {0}.".format(
tz.__getitem__(index=0).to_TOA_list()[0].obs
)
assert "clkcorr" not in tz.__getitem__(index=0).table["flags"][0], failstr


def test_tempo2_consistency(toas, model):
"""Test if derived BATs and absolute phases are consisent with those of tempo2."""
bary_pint = model.get_barycentric_toas(toas)
phase_pint = model.phase(toas, abs_phase=True)[1]
phase_pint = phase_pint % 1

# Tolerance of 1us, it is smaller than the typical 27us
# from TT(BIPM).
atol = (1.0 * u.us).to(u.d)

bary_test = u.isclose(bary_tempo2, bary_pint, rtol=0, atol=atol)

# The same for phases, times the frequency in the sample
# par file.
atol = (1.0 * u.us) * (100 * u.Hz)

phase_test = u.isclose(phase_tempo2, phase_pint, rtol=0, atol=atol)

delta_bary = bary_tempo2 - bary_pint
delta_phase = phase_tempo2 - phase_pint

for db, dp in zip(delta_bary, delta_phase):
dbnano = db.to(u.ns)
print(f"BAT delta: {db:+.6e}, {dbnano:+.6e}, phase {dp:+.6e}")

for db, dp, BARY_TEST, PHASE_TEST in zip(
delta_bary, delta_phase, bary_test, phase_test
):
assert BARY_TEST, f"BAT delta is beyond 1us: {db:+.6e}"
assert PHASE_TEST, f"PHASE delta is beyons 1us: {dp:+.6e}"
Loading