|
| 1 | +"""System and frequency dependent delays to model profile evolution.""" |
| 2 | + |
| 3 | +import re |
| 4 | +from warnings import warn |
| 5 | + |
| 6 | +import astropy.units as u |
| 7 | +import numpy as np |
| 8 | + |
| 9 | +from pint.models.parameter import boolParameter, maskParameter |
| 10 | +from pint.models.timing_model import DelayComponent |
| 11 | + |
| 12 | +fdjump_max_index = 20 |
| 13 | + |
| 14 | + |
| 15 | +class FDJump(DelayComponent): |
| 16 | + """A timing model for system-dependent frequency evolution of pulsar |
| 17 | + profiles. |
| 18 | +
|
| 19 | + This model expresses the delay as a polynomial function of the |
| 20 | + observing frequency/logarithm of observing frequency in the SSB frame. |
| 21 | + This is intended to compensate for the delays introduced by frequency-dependent |
| 22 | + profile structure when a different profiles are used for different systems. |
| 23 | +
|
| 24 | + The default behavior is to have FDJUMPs as polynomials of the observing |
| 25 | + frequency (rather than log-frequency). This is different from the convention |
| 26 | + used for global FD parameters. This choice is made to be compatible with tempo2. |
| 27 | + This is controlled using the FDJUMPLOG parameter. "FDJUMPLOG Y" may not be |
| 28 | + tempo2-compatible. |
| 29 | +
|
| 30 | + Note |
| 31 | + ---- |
| 32 | + FDJUMPs have two indices: the polynomial/FD/prefix index and the system/mask |
| 33 | + index. i.e., they have properties of both maskParameters such as JUMPs and |
| 34 | + prefixParameters such as FDs. There is currently no elegant way in PINT to implement |
| 35 | + such parameters due to the way parameter indexing is implemented; there is no way to |
| 36 | + distinguish between mask and prefix indices. |
| 37 | +
|
| 38 | + Hence, they are implemented here as maskParameters as a stopgap measure. |
| 39 | + This means that there must be an upper limit for the FD indices. This is controlled |
| 40 | + using the `pint.models.fdjump.fdjump_max_index` variable, and is 20 by default. |
| 41 | + Note that this is strictly a limitation of the implementation and not a property |
| 42 | + of FDJUMPs themselves. |
| 43 | +
|
| 44 | + FDJUMPs appear in tempo2-format par files as "FDJUMPp", where p is the FD index. |
| 45 | + The mask index is not explicitly mentioned in par files similar to JUMPs. |
| 46 | + PINT understands both "FDJUMPp" and "FDpJUMP" as the same parameter in par files, |
| 47 | + but the internal representation is always "FDpJUMPq", where q is the mask index. |
| 48 | +
|
| 49 | + PINT understands 'q' as the mask parameter just fine, but the identification of 'p' |
| 50 | + as the prefix parameter is done in a hacky way. |
| 51 | +
|
| 52 | + This implementation may be overhauled in the future. |
| 53 | +
|
| 54 | + Parameters supported: |
| 55 | +
|
| 56 | + .. paramtable:: |
| 57 | + :class: pint.models.fdjump.FDJump |
| 58 | + """ |
| 59 | + |
| 60 | + register = True |
| 61 | + category = "fdjump" |
| 62 | + |
| 63 | + def __init__(self): |
| 64 | + super().__init__() |
| 65 | + |
| 66 | + # Matches "FDpJUMPq" where p and q are integers. |
| 67 | + self.param_regex = re.compile("^FD(\\d+)JUMP(\\d+)") |
| 68 | + |
| 69 | + self.add_param( |
| 70 | + boolParameter( |
| 71 | + name="FDJUMPLOG", |
| 72 | + value=False, |
| 73 | + description="Whether to use log-frequency (Y) or linear-frequency (N) for computing FDJUMPs.", |
| 74 | + ) |
| 75 | + ) |
| 76 | + for j in range(1, fdjump_max_index + 1): |
| 77 | + self.add_param( |
| 78 | + maskParameter( |
| 79 | + name=f"FD{j}JUMP", |
| 80 | + units="second", |
| 81 | + description=f"System-dependent FD parameter of polynomial index {j}", |
| 82 | + ) |
| 83 | + ) |
| 84 | + |
| 85 | + self.delay_funcs_component += [self.fdjump_delay] |
| 86 | + |
| 87 | + def setup(self): |
| 88 | + super().setup() |
| 89 | + |
| 90 | + self.fdjumps = [ |
| 91 | + mask_par |
| 92 | + for mask_par in self.get_params_of_type("maskParameter") |
| 93 | + if self.param_regex.match(mask_par) |
| 94 | + ] |
| 95 | + |
| 96 | + for fdj in self.fdjumps: |
| 97 | + # prevents duplicates from being added to phase_deriv_funcs |
| 98 | + if fdj in self.deriv_funcs.keys(): |
| 99 | + del self.deriv_funcs[fdj] |
| 100 | + self.register_deriv_funcs(self.d_delay_d_FDJUMP, fdj) |
| 101 | + |
| 102 | + def get_fd_index(self, par): |
| 103 | + """Extract the FD index from an FDJUMP parameter name. In a parameter name |
| 104 | + "FDpJUMPq", p is the FD/prefix index and q is the mask index. |
| 105 | +
|
| 106 | + Parameters |
| 107 | + ---------- |
| 108 | + par: Parameter name (str) |
| 109 | +
|
| 110 | + Returns |
| 111 | + ------- |
| 112 | + FD index (int) |
| 113 | + """ |
| 114 | + if m := self.param_regex.match(par): |
| 115 | + return int(m.groups()[0]) |
| 116 | + else: |
| 117 | + raise ValueError( |
| 118 | + f"The given parameter {par} does not correspond to an FDJUMP." |
| 119 | + ) |
| 120 | + |
| 121 | + def get_freq_y(self, toas): |
| 122 | + """Get frequency or log-frequency in GHz based on the FDJUMPLOG value. |
| 123 | + Returns (freq/1_GHz) if FDJUMPLOG==N and log(freq/1_GHz) if FDJUMPLOG==Y. |
| 124 | + Any non-finite values are replaced by zero. |
| 125 | +
|
| 126 | + Parameters |
| 127 | + ---------- |
| 128 | + toas: pint.toa.TOAs |
| 129 | +
|
| 130 | + Returns |
| 131 | + ------- |
| 132 | + (freq/1_GHz) or log(freq/1_GHz) depending on the value of FDJUMPLOG (float). |
| 133 | + """ |
| 134 | + tbl = toas.table |
| 135 | + try: |
| 136 | + freq = self._parent.barycentric_radio_freq(toas) |
| 137 | + except AttributeError: |
| 138 | + warn("Using topocentric frequency for frequency dependent delay!") |
| 139 | + freq = tbl["freq"] |
| 140 | + |
| 141 | + y = ( |
| 142 | + np.log(freq.to(u.GHz).value) |
| 143 | + if self.FDJUMPLOG.value |
| 144 | + else freq.to(u.GHz).value |
| 145 | + ) |
| 146 | + non_finite = np.invert(np.isfinite(y)) |
| 147 | + y[non_finite] = 0.0 |
| 148 | + |
| 149 | + return y |
| 150 | + |
| 151 | + def fdjump_delay(self, toas, acc_delay=None): |
| 152 | + """Calculate frequency dependent delay. |
| 153 | +
|
| 154 | + If FDJUMPLOG is Y, use the following expression (similar to global FD parameters): |
| 155 | +
|
| 156 | + FDJUMP_delay = sum_i(c_i * (log(obs_freq/1GHz))^i) |
| 157 | +
|
| 158 | + If FDJUMPLOG is N, use the following expression (same as in tempo2, default): |
| 159 | +
|
| 160 | + FDJUMP_delay = sum_i(c_i * (obs_freq/1GHz)^i) |
| 161 | + """ |
| 162 | + y = self.get_freq_y(toas) |
| 163 | + |
| 164 | + delay = np.zeros_like(y) |
| 165 | + for fdjump in self.fdjumps: |
| 166 | + fdj = getattr(self, fdjump) |
| 167 | + if fdj.quantity is not None: |
| 168 | + mask = fdj.select_toa_mask(toas) |
| 169 | + ymask = y[mask] |
| 170 | + fdidx = self.get_fd_index(fdjump) |
| 171 | + fdcoeff = fdj.value |
| 172 | + delay[mask] += fdcoeff * ymask**fdidx |
| 173 | + |
| 174 | + return delay * u.s |
| 175 | + |
| 176 | + def d_delay_d_FDJUMP(self, toas, param, acc_delay=None): |
| 177 | + """Derivative of delay w.r.t. FDJUMP parameters.""" |
| 178 | + assert ( |
| 179 | + bool(self.param_regex.match(param)) |
| 180 | + and hasattr(self, param) |
| 181 | + and getattr(self, param).quantity is not None |
| 182 | + ), f"{param} is not present in the FDJUMP model." |
| 183 | + |
| 184 | + y = self.get_freq_y(toas) |
| 185 | + mask = getattr(self, param).select_toa_mask(toas) |
| 186 | + ymask = y[mask] |
| 187 | + fdidx = self.get_fd_index(param) |
| 188 | + |
| 189 | + delay_derivative = np.zeros_like(y) |
| 190 | + delay_derivative[mask] = ymask**fdidx |
| 191 | + |
| 192 | + return delay_derivative * u.dimensionless_unscaled |
| 193 | + |
| 194 | + def print_par(self, format="pint"): |
| 195 | + par = super().print_par(format) |
| 196 | + |
| 197 | + if format != "tempo2": |
| 198 | + return par |
| 199 | + |
| 200 | + for fdjump in self.fdjumps: |
| 201 | + if getattr(self, fdjump).quantity is not None: |
| 202 | + j = self.get_fd_index(fdjump) |
| 203 | + par = par.replace(f"FD{j}JUMP", f"FDJUMP{j}") |
| 204 | + |
| 205 | + return par |
0 commit comments