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

PDA Updater #913

Merged
merged 7 commits into from
May 13, 2024
Merged

PDA Updater #913

merged 7 commits into from
May 13, 2024

Conversation

jmbarr
Copy link
Contributor

@jmbarr jmbarr commented Dec 20, 2023

Provides an updater based on probabilistic data association. See the documentation in the Updaters section, or have a look at this gist: https://gist.github.com/jmbarr/92dc83e28c04026136d4f8706a1159c1

@jmbarr jmbarr requested a review from a team as a code owner December 20, 2023 09:06
@jmbarr jmbarr requested review from nperree-dstl and spike-dstl and removed request for a team December 20, 2023 09:06
Comment on lines +187 to +190
posterior_mean += kalman_gain @ sum_of_innovations
posterior_covariance += \
kalman_gain @ (sum_of_weighted_cov - sum_of_innovations @ sum_of_innovations.T) \
@ kalman_gain.T
Copy link
Contributor

Choose a reason for hiding this comment

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

Do these lines work if line 162 is False? I can't work out where posterior_mean, posterior_covariance or kalman_gain would come from in that case.

Copy link
Contributor

Choose a reason for hiding this comment

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

Ignore me, I think I've got it now - the first hypothesis always a MissedDetection hypothesis right?

Copy link
Contributor Author

Choose a reason for hiding this comment

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

the first hypothesis always a MissedDetection hypothesis right?

It seems to be but that's a matter for the associator.associate() method. I think this should work whatever the order; it exploits the fact that not hypothesis returns True for a missed detection (not sure why that is) and calculates the required mean and covariance, and attaches a 0-innovation.

Comment on lines +177 to +182
# probably exists a less clunky way of doing this using exists() or overwritten +=
# All these floats should be redundant if/when the bug in Probability.__mult__() is
# fixed.
if n == 0:
sum_of_innovations = float(hypothesis.probability) * innovation
sum_of_weighted_cov = float(hypothesis.probability) * (innovation @ innovation.T)
Copy link
Member

Choose a reason for hiding this comment

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

Set these two sum variables outside the loop to 0?

Copy link
Member

@sdhiscocks sdhiscocks Jan 3, 2024

Choose a reason for hiding this comment

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

Wasn't too sure if the proposed would work, but doing a simple test, looks like it should.

>>> a = 0
>>> a += np.array([[-5], [2]])
>>> a
array([[-5],
       [ 2]])

Copy link
Contributor Author

@jmbarr jmbarr Jan 3, 2024

Choose a reason for hiding this comment

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

So,
sum_of_innovations = 0
sum_of_innovations += hypothesis.probability * innovation
?

Copy link
Contributor Author

Choose a reason for hiding this comment

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

So, sum_of_innovations = 0 sum_of_innovations += hypothesis.probability * innovation ?

But won't matter if :class:Probability is corrected.

Copy link
Member

Choose a reason for hiding this comment

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

diff --git a/stonesoup/updater/probability.py b/stonesoup/updater/probability.py
index c471e514..c6495111 100644
--- a/stonesoup/updater/probability.py
+++ b/stonesoup/updater/probability.py
@@ -156,6 +156,8 @@ class PDAUpdater(ExtendedKalmanUpdater):
             The covariance of the reduced/single Gaussian
         """
 
+        sum_of_innovations = 0
+        sum_of_weighted_cov = 0
         for n, hypothesis in enumerate(hypotheses):
             # Check for the existence of an associated measurement. Because of the way the
             # hypothesis is constructed, you can do this:
@@ -174,15 +176,10 @@ class PDAUpdater(ExtendedKalmanUpdater):
                 innovation = hypothesis.measurement.state_vector - \
                              hypothesis.measurement_prediction.state_vector
 
-            # probably exists a less clunky way of doing this using exists() or overwritten +=
             # All these floats should be redundant if/when the bug in Probability.__mult__() is
             # fixed.
-            if n == 0:
-                sum_of_innovations = float(hypothesis.probability) * innovation
-                sum_of_weighted_cov = float(hypothesis.probability) * (innovation @ innovation.T)
-            else:
-                sum_of_innovations += float(hypothesis.probability) * innovation
-                sum_of_weighted_cov += float(hypothesis.probability) * (innovation @ innovation.T)
+            sum_of_innovations += float(hypothesis.probability) * innovation
+            sum_of_weighted_cov += float(hypothesis.probability) * (innovation @ innovation.T)
 
         posterior_mean += kalman_gain @ sum_of_innovations
         posterior_covariance += \

@sglvladi
Copy link
Collaborator

Looking at this made me think whether we could make this more modular by using a KalmanUpdater (or subclass) via composition, rather than inheritance:

  • PDAUpdater inherits from Updater and accepts a KalmanUpdater as an updater property, which it then uses to either update() or predict_measurement(), as per lines 127 and 164
  • The class could also be renamed to KalmanPDAUpdater/GaussianPDAUpdater to reflect the fact that it applies specifically to Kalman/Gaussian cases.

The above would pave the way for a ParticlePDAUpdater/ELPFUpdater that could be structured similarly to perform the PDA update for particle filters.

@sdhiscocks
Copy link
Member

Looking at this made me think whether we could make this more modular by using a KalmanUpdater (or subclass) via composition, rather than inheritance:

  • PDAUpdater inherits from Updater and accepts a KalmanUpdater as an updater property, which it then uses to either update() or predict_measurement(), as per lines 127 and 164
  • The class could also be renamed to KalmanPDAUpdater/GaussianPDAUpdater to reflect the fact that it applies specifically to Kalman/Gaussian cases.

The above would pave the way for a ParticlePDAUpdater/ELPFUpdater that could be structured similarly to perform the PDA update for particle filters.

Something similar was done for the DIEKF in #891: https://github.com/dstl/Stone-Soup/blob/e27fc3639d791d2441b6c4b8572c57932651deb3/stonesoup/updater/iterated.py

@sdhiscocks sdhiscocks merged commit 0089d29 into main May 13, 2024
6 checks passed
@sdhiscocks sdhiscocks deleted the PDAUpdater2 branch May 13, 2024 08:58
@sdhiscocks
Copy link
Member

Per @sglvladi comment, I think it would be good to change structure to support more flexible structure, but I'll merge this for now.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Projects
None yet
Development

Successfully merging this pull request may close these issues.

5 participants