Skip to content

Commit 6dfe23a

Browse files
committed
added a little to the readme and references to issue 178 in comments
'
1 parent 0891a55 commit 6dfe23a

2 files changed

Lines changed: 9 additions & 7 deletions

File tree

README.md

Lines changed: 5 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -24,7 +24,7 @@ Python methods for numerical differentiation of noisy data, including multi-obje
2424

2525
## Introduction
2626

27-
PyNumDiff is a Python package that implements various methods for computing numerical derivatives of noisy data, which can be a critical step in developing dynamic models or designing control. There are seven different families of methods implemented in this repository:
27+
PyNumDiff is a Python package that implements many methods for computing numerical derivatives and smooth estimates of noisy data, which can be a critical step in developing dynamic models or designing control. There are seven different families of methods implemented in this repository:
2828

2929
1. prefiltering followed by finite difference calculation
3030
2. iterated finite differencing
@@ -34,9 +34,11 @@ PyNumDiff is a Python package that implements various methods for computing nume
3434
6. generalized Kalman smoothing
3535
7. local approximation with linear model
3636

37-
For a full list, explore modules in the [Sphinx documentation](https://pynumdiff.readthedocs.io/master/), or read section 7 of our [Taxonomy Paper](https://arxiv.org/abs/2512.09090).
37+
All are ultimately smoothing with similar runtime and accuracy, but some have situational advantages over others: For example, `robustdiff` is specialized to handle outliers; `splinediff`, `polydiff`, `rtsdiff`, and `robustdiff` can handle missing data; `splinediff`, `polydiff`, `rbfdiff`, `rtsdiff`, and `robustdiff` can handle irregularly-spaced data; and `rtsdiff` can handle inputs on a wrapped domain, like angles. All methods can accept blocks of multidimensional data, differentiating all vectors along the dimension given by the `axis` parameter.
3838

39-
Most of these methods have multiple parameters, so we take a principled approach and propose a multi-objective optimization framework for choosing parameters that minimize a loss function to balance the faithfulness and smoothness of the derivative estimate. For more details, refer to [this paper](https://doi.org/10.1109/ACCESS.2020.3034077).
39+
For a full list and comparison, see section 7 of our [Taxonomy Paper](https://arxiv.org/abs/2512.09090) and explore modules in the [Sphinx documentation](https://pynumdiff.readthedocs.io/master/).
40+
41+
All methods have hyperparameters, so we take a principled approach and propose a multi-objective optimization framework for choosing settings that minimize a loss function to balance the faithfulness and smoothness of the derivative estimate. For more details, refer to [this paper](https://doi.org/10.1109/ACCESS.2020.3034077).
4042

4143
## Installing
4244

pynumdiff/kalman_smooth.py

Lines changed: 4 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -25,8 +25,8 @@ def kalman_filter(y, xhat0, P0, A, Q, C, R, B=None, u=None, save_P=True, innovat
2525
:param bool save_P: whether to save history of error covariance and a priori state estimates, used with rts
2626
smoothing but nonstandard to compute for ordinary filtering
2727
:param callable innovation_fn: optional function taking measurements and predicted measurements and returning the innovation.
28-
When :code:`None`, traditional subtraction is used. This is exposed to handle cases like wrapped domains, where alternative
29-
displacement measures may be more appropriate. See e.g. the function passed by :code:`rtsdiff` with :code:`circular=True`.
28+
When :code:`None`, traditional subtraction is used. This is primarily exposed to handle angles, which have a wrapped domain,
29+
so alternative displacement measure :code:`lambda y, pred: (y - pred + np.pi) % (2*np.pi) - np.pi` is more appropriate.
3030
3131
:return: - **xhat_pre** (np.array) -- a priori estimates of xhat, with axis=0 the batch dimension, so xhat[n] gets the nth step
3232
- **xhat_post** (np.array) -- a posteriori estimates of xhat
@@ -147,7 +147,7 @@ def rtsdiff(x, dt_or_t, order, log_qr_ratio, forwardbackward=False, axis=0, circ
147147
Q_d[n] = eM[:order+1, order+1:] @ A_d[n].T
148148
if forwardbackward: A_d_bwd = np.linalg.inv(A_d[::-1]) # properly broadcasts, taking inv of each stacked 2D array
149149

150-
innovation_fn = (lambda y, pred: (y - pred + np.pi) % (2*np.pi) - np.pi) if circular else None # wrap innovation to [-pi, pi] for circular variables
150+
innovation_fn = (lambda y, pred: (y - pred + np.pi) % (2*np.pi) - np.pi) if circular else None # optionall wrap innovation to [-pi, pi], see #178
151151

152152
x_hat = np.empty_like(x); dxdt_hat = np.empty_like(x)
153153
if forwardbackward: w = np.linspace(0, 1, N) # weights used to combine forward and backward results
@@ -171,7 +171,7 @@ def rtsdiff(x, dt_or_t, order, log_qr_ratio, forwardbackward=False, axis=0, circ
171171
x_hat[s] = x_hat[s] * w + xhat_smooth[:, 0][::-1] * (1-w)
172172
dxdt_hat[s] = dxdt_hat[s] * w + xhat_smooth[:, 1][::-1] * (1-w)
173173

174-
if circular: x_hat = (x_hat + np.pi) % (2*np.pi) - np.pi # wrap output to match the input domain
174+
if circular: x_hat = (x_hat + np.pi) % (2*np.pi) - np.pi # wrap output to match the input domain, see #178
175175
return x_hat, dxdt_hat
176176

177177

0 commit comments

Comments
 (0)