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

Add method to create events with interacting charge clouds #460

Merged
merged 1 commit into from
Feb 13, 2025

Conversation

fhagemann
Copy link
Collaborator

There are several ways of creating a multi-site Event, depending on whether the individual charge clouds should interact amongst each other or not.

Example:
Event([[p1], [p2,p3], [p4,p5]], energies) --> in this case, the charge cloud at p1 would drift independently, the charge clouds at p2 and p3 could interact with each other and the ones at p4 and p5 could interact with each other.
There is a convenience function that transforms Event([p1, p2]) to Event([[p1], [p2]]), meaning that p1 and p2 cannot interact with each other. If they should, the syntax would have to be Event([[p1, p2]]) (note the extra pair of brackets).

Visual example of what the extra set of brackets would cause:

using SolidStateDetectors
using Unitful
using Plots
T = Float32

sim = Simulation{T}(SSD_examples[:InvertedCoax])
simulate!(sim)

evt1 = Event([CartesianPoint{T}(0.02,0,0.05), CartesianPoint{T}(0.02,0.0,0.051)], [1u"GeV", 1u"GeV"])
simulate!(evt1, sim, Δt = 1u"ns", max_nsteps = 5000, self_repulsion = true)
p1 = plot(sim.detector)
plot!(p1, evt1.drift_paths)

evt2 = Event([[CartesianPoint{T}(0.02,0,0.05), CartesianPoint{T}(0.02,0.0,0.051)]], [[1u"GeV", 1u"GeV"]])
simulate!(evt2, sim, Δt = 1u"ns", max_nsteps = 5000, self_repulsion = true)
p2 = plot(sim.detector)
plot!(p2, evt2.drift_paths)
    
plot(p1, p2, size = (1000,500), legend = false)

image

Left: the charge clouds do not repell each other. Right: the charge clouds "see" and repell each other.

When wanting to create charge cloud, a third argument N (number of charges) can be passed to Event.
This works with the syntax in evt1 (independent charge clouds), but not with the syntax in evt2 (combined charge clouds).

I added support for the latter syntax in this PR:

using SolidStateDetectors
using Unitful
using Plots
T = Float32

sim = Simulation{T}(SSD_examples[:InvertedCoax])
simulate!(sim)

evt1 = Event([CartesianPoint{T}(0.02,0,0.05), CartesianPoint{T}(0.02,0.0,0.051)], [1u"GeV", 1u"GeV"], 5)
simulate!(evt1, sim, Δt = 1u"ns", max_nsteps = 5000, self_repulsion = true)
p1 = plot(sim.detector)
plot!(p1, evt1.drift_paths)

evt2 = Event([[CartesianPoint{T}(0.02,0,0.05), CartesianPoint{T}(0.02,0.0,0.051)]], [[1u"GeV", 1u"GeV"]], 5)
simulate!(evt2, sim, Δt = 1u"ns", max_nsteps = 5000, self_repulsion = true)
p2 = plot(sim.detector)
plot!(p2, evt2.drift_paths)
    
plot(p1, p2, size = (1000,500), legend = false)

image

With the new method, if one were to split one deposit into two, the resulting waveforms would look the same (whereas with the independent charge clouds, effects of self_repulsion would be underestimated):

using SolidStateDetectors
using Unitful
using Plots
T = Float32

sim = Simulation{T}(SSD_examples[:InvertedCoax])
simulate!(sim)

pos = CartesianPoint{T}(0.02,0,0.05)
E = T(1)u"GeV"

evt0 = Event([pos], [E], 50)
simulate!(evt0, sim, Δt = 1u"ns", max_nsteps = 5000, self_repulsion = true)

evt1 = Event([pos, pos], [E/2, E/2], 50)
simulate!(evt1, sim, Δt = 1u"ns", max_nsteps = 5000, self_repulsion = true)

evt2 = Event([[pos, pos]], [[E/2, E/2]], 50)
simulate!(evt2, sim, Δt = 1u"ns", max_nsteps = 5000, self_repulsion = true)

plot(evt0.waveforms[1],  label = "One hit")
plot!(evt1.waveforms[1], label = "Two hits, independent")
plot!(evt2.waveforms[1], label = "Two hits, interacting")

image

@fhagemann fhagemann added the enhancement Improvement of existing features label Feb 13, 2025
@hervasa2
Copy link
Collaborator

Amazing! Works as expected.

@fhagemann fhagemann merged commit 63e61f4 into main Feb 13, 2025
9 checks passed
@fhagemann fhagemann deleted the charge_cloud_events branch February 13, 2025 22:39
@tdixon97
Copy link

This is interesting, but if we had many charge clouds interacting would the simulation get unreasonably slow?

@fhagemann
Copy link
Collaborator Author

That depends on the number of charges you set. The current implementation scales with O(N^2).

@tdixon97
Copy link

tdixon97 commented Feb 18, 2025

But is N the total number across all the clouds ? I guess so. Probably we need to consider treating some clouds as interacting (those close) and some as non-interacting? Is that currently possible?

@fhagemann
Copy link
Collaborator Author

fhagemann commented Feb 18, 2025

Right now, N is applied to every charge cloud, not for all the clouds.
With this PR, you can cluster the initial points passed to Event.

If you want p1 and p2 to interact, but not p3, you can MANUALLY call
Event([[p1, p2], [p3]], [[E1, E2], E3]]). But we (@apmypb) are currently working on providing a grouping algorithm that would do the grouping AUTOMATICALLY based on some given max_distance above which the charge clouds should not interact.

@tdixon97
Copy link

Ok that sounds quite interesting for applications with a lot of steps (gammas)

@fhagemann
Copy link
Collaborator Author

Ok that sounds quite interesting for applications with a lot of steps (gammas)

@tdixon97 @hervasa2 Please have a look at #461, which should allow for grouping hits automatically based on a given group_distance

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
enhancement Improvement of existing features
Projects
None yet
Development

Successfully merging this pull request may close these issues.

3 participants