Skip to content

Add time rescaling to inhomogeneous#81

Merged
rsenne merged 6 commits into
mainfrom
time_change_inhomogeneous_
May 10, 2026
Merged

Add time rescaling to inhomogeneous#81
rsenne merged 6 commits into
mainfrom
time_change_inhomogeneous_

Conversation

@rsenne
Copy link
Copy Markdown
Collaborator

@rsenne rsenne commented May 9, 2026

Resolves #78

  • Add time change method for inhomogeneous pps
  • Tests time change method
  • Tests for hypothesis tests and inhomogeneous
  • Adds to inhomogeneous docs

@codecov
Copy link
Copy Markdown

codecov Bot commented May 9, 2026

Codecov Report

✅ All modified and coverable lines are covered by tests.

📢 Thoughts on this report? Let us know!

Copy link
Copy Markdown

Copilot AI left a comment

Choose a reason for hiding this comment

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

Pull request overview

Adds a time-rescaling (“time change”) transform for InhomogeneousPoissonProcess so existing goodness-of-fit tests (KS vs Uniform / Exponential) can be applied after compensator-based rescaling, with accompanying tests and documentation updates.

Changes:

  • Implement time_change(h, pp::InhomogeneousPoissonProcess) using the integrated intensity (compensator).
  • Add unit tests for time rescaling on polynomial and custom intensities, plus integration into hypothesis-test workflows.
  • Extend the inhomogeneous Poisson documentation example with a Monte Carlo KS goodness-of-fit section.

Reviewed changes

Copilot reviewed 4 out of 4 changed files in this pull request and generated 3 comments.

File Description
src/univariate/poisson/inhomogeneous/inhomogeneous_poisson_process.jl Adds time_change for inhomogeneous Poisson processes via integrated_intensity.
test/inhomogeneous_poisson_process.jl Adds direct tests for correctness of the time-change transform (incl. non-zero tmin, empty history, custom intensity).
test/hypothesis_tests.jl Adds tests exercising KSDistance/MonteCarloTest integration for inhomogeneous Poisson via time_change.
docs/examples/Inhomogeneous.jl Documents goodness-of-fit testing via time rescaling and Monte Carlo KS tests.

💡 Add Copilot custom instructions for smarter, more guided reviews. Learn how to get started.

Comment thread test/hypothesis_tests.jl Outdated
Comment thread test/inhomogeneous_poisson_process.jl Outdated
rsenne and others added 2 commits May 8, 2026 23:09
…cess.jl

Co-authored-by: github-actions[bot] <41898282+github-actions[bot]@users.noreply.github.com>
Comment on lines +125 to +139
new_tmax = integrated_intensity(f, h.tmin, h.tmax, config)
T = typeof(new_tmax)
new_times = T[integrated_intensity(f, h.tmin, t, config) for t in h.times]
# `new_tmax` and each `new_times[i]` come from independent quadrature calls,
# so solver / floating-point error can let a transformed time slightly exceed
# `new_tmax`. Widen `new_tmax` to the observed maximum and skip History's
# bounds check so events are never silently dropped.
if !isempty(new_times)
new_tmax = max(new_tmax, last(new_times) * (one(T) + eps(T)))
end
return History(;
times=new_times, tmin=zero(T), tmax=new_tmax, marks=h.marks, check_args=false
)
end

Copy link
Copy Markdown
Owner

Choose a reason for hiding this comment

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

If we have to check if the order was changed for new_tmax, then isn't it possible that the order of any pair of events got switched? In this case, should the History be initialized with check_args set to true?

Copy link
Copy Markdown
Collaborator Author

Choose a reason for hiding this comment

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

You are right. My mistake. After reviewing the code of History i think this speaks to an underlying issue with our logic for input handling. I think History might be better of as a validator and not a fixer. What i mean by this is I feel history might be better of only qualifying if the history times meet the proper checks (i.e., bounded properly and sorted).

My thinking for this is that there are three reasonable cases unsorted data happens: i.) a user has unsorted data and needs to sort it ii.) a user has an actual error in their data they don't know about iii) we have numerical issues that are introduced from independent quadrature calls.

Right now the behavior of History basically conflates all three as equally valid/equally reasonable occurrences. I don't think we should try to reason about what the case of the above is except when we know for sure what is happening like in numerical issue introduced by quadrature.

On our side my thought is we handle these numerical issues and then use History as a way to validate our own code works as expected (e.g., by sorting events and also expanding the tmax). In the case of the user i think we should just throw an error. That way they have to opt-in to any changes of their data and or reciprocally opt-out.

This is getting a bit off-topic for this specific case, but I think this is something that might be at the least thinking about more concretely. We could easily supply some easy helpers so that a user knows they are opting into these sort of data fixes e.g., sort_history! and restrict. This would be a breaking change but we aren't at major version 1 yet so semantic versioning doesn't require us to notify users but we could at least add a deprecation cycle if we wanted to.

Hope this makes sense--let me know if it doesn't! I can also move this over to a new issue.

Copy link
Copy Markdown
Owner

Choose a reason for hiding this comment

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

My answer to the History stuff will be in the main thread. But about the ordering, since we will need to call quadrature n+1 times anyway, can't we accumulate the integral of interarrival times instead of calculating the integral of each event time? I would guess quadrature would not return negative results for non negative functions, therefore this would not change the ordering.

@JoseKling
Copy link
Copy Markdown
Owner

Great. I was planning on implementing this next.
Just check if it makes sense what I commented, but all good.

@rsenne rsenne requested a review from JoseKling May 9, 2026 14:16
@rsenne
Copy link
Copy Markdown
Collaborator Author

rsenne commented May 9, 2026

I re-requested review as I made further changes. I edited the behavior per your comment (see the response I left i think it's worth thinking about), and also made JET runtime conditional as the constant failure on pre-release was annoying me lol. This PR got a little noisy as a result so I'd probably squash merge it

@JoseKling
Copy link
Copy Markdown
Owner

About History being a validator or a fixer, I thought about it a while.
I think this validator idea makes sense, but I am not sure we should use it in every case. If we are performing a lot of simulations, like for goodness-of-fit tests, we don't want to check the ordering in every simulation, the algorithm should already guarantee that it is correct.
So what I was thinking is that we could keep History as it is, but throwing a warning if the event times are not ordered. This would validate the data, but fix it at the same time.
Then internally we should always use check_args=False. The methods themselves should guarantee compatibility.
I guess we loose some safety by not using the validation, but ideally the algorithms obviously produce valid entries or have been well tested.
Not 100% on this, so let me know what you think.

Comment thread test/runtests.jl
Comment on lines +23 to +31
# Skip JET on Julia pre-releases (where JET typically hasn't caught up
# yet and there is no compatible JET version to resolve against). JET is
# not listed in test/Project.toml's [deps] for the same reaso, having
# it there would make `Pkg.test` fail at the resolution step on
# prereleases, before this guard ever runs. Install on demand only when
# we're on a stable Julia.
if isempty(VERSION.prerelease)
Pkg.add("JET")
@eval using JET
Copy link
Copy Markdown
Owner

Choose a reason for hiding this comment

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

great. That was annoying ; )

@rsenne
Copy link
Copy Markdown
Collaborator Author

rsenne commented May 10, 2026

On the cumulative inter-arrival idea: i think this is exactly right and we should do this in another PR. If we integrate over [tmin, t_1], [t_1, t_2], ..., [t_n, tmax] and cumsum it, ordering is monotone by construction and new_tmax is literally last(new_times) + last(deltas), so its automatically >= last(new_times). Same number of quadrature calls, but we get rid of both the order swap risk and the widening hack in one go. One paranoid caveat is that adaptive quadrature can in theory return a tiny negative result on a non negative integrand in pathological cases, so we could just clamp with max.(zero(T), deltas) before the cumsum. Costs nothing and should be bulletproof.

On keeping History as a fixer with a warning, i see the appeal, especially the worry about MC overhead, but im not certain this would cost us much. The actual validation is just issorted plus a searchsortedfirst, which is nanoseconds next to even one quadrature call so im skeptical it would show up in a profile. And @warn inside a MC loop has its own problem, either it fires every iteration and floods the logs (and @warn isnt free either), or it doesnt fire and we're back to silently fixing things. The thing that worries me most is that "internal callers use check_args=false and trust the algorithm" is the pattern that bit us in time_change.

Nonetheless, I'm going to merge this PR and I can start a discussion on these points as i view them as separate issues whether added or not

@rsenne rsenne merged commit 7124e07 into main May 10, 2026
6 checks passed
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

None yet

Projects

None yet

Development

Successfully merging this pull request may close these issues.

Add time_change method to inhomogeneous so we can use time rescaling

3 participants