fix: align post-boundary optical timing with Geant4 11.4+#260
Conversation
779b7be to
cff78e4
Compare
Geant4 11.4 substantially refactored G4OpBoundaryProcess, including the dielectric-dielectric boundary handling path. In our simg4ox and GPURaytrace workflows, Geant4 11.4+ is the first release family where timing-sensitive optical comparisons diverge from the Geant4 11.3 baselines. That makes post-boundary timing correctness explicit in a way our previous implementation did not handle. This branch fixes two gaps in our handling of post-boundary optical timing: 1. The Opticks GPU path advanced the next segment using the pre-boundary medium group velocity even after transmission into material2. 2. The lightweight Geant4-side tracking/recording path used by simg4ox did not explicitly enable UseGivenVelocity(true) for optical photons, so the track was not forced to follow the boundary-selected velocity during post-boundary recording. Fix the GPU side by carrying the active-medium timing state through qsim: - store both material1 and material2 group velocities in sstate - track current_material_index/current_group_velocity in sctx - use current_group_velocity in propagate_to_boundary - switch the carried material/velocity in propagate_at_boundary after reflection/transmission - refresh the carried velocity after each boundary lookup when the boundary orientation still matches the active medium - initialize the new fields in the QSim test/debug kernels Fix the Geant4-side CPU recording path by enabling UseGivenVelocity(true) for optical photons in g4app tracking action, so the recorded post-boundary timing follows the velocity selected by Geant4. Because hit merging in this codebase is time-bucket based, correcting post-boundary timing changes timing-sensitive observables and can change merged hit counts. The follow-up test updates on this branch therefore make expected residual diffs and hit-count baselines depend on the Geant4 release family instead of treating the new values as regressions.
cff78e4 to
c909f7d
Compare
c909f7d to
7b14dea
Compare
|
@codex review |
There was a problem hiding this comment.
Pull request overview
Aligns Opticks GPU and Geant4-side CPU recording with Geant4 11.4+ optical boundary timing behavior by carrying/using the post-boundary (active medium) group velocity across boundary interactions, and updating tests to use Geant4-series-specific baselines.
Changes:
- GPU: extend
sstate/sctxand updateqsimpropagation to use a carried “current” group velocity and switch it correctly on reflection/transmission. - CPU recording: force optical tracks to honor boundary-proposed velocity via
UseGivenVelocity(true)in tracking action. - Tests: make hit-count and record-diff expectations depend on Geant4 release series (
11.3vs11.4+).
Reviewed changes
Copilot reviewed 9 out of 9 changed files in this pull request and generated 2 comments.
Show a summary per file
| File | Description |
|---|---|
| tests/test_GPURaytrace.sh | Selects expected hit-count baselines by detected Geant4 series. |
| tests/compare_ab.py | Uses detected Geant4 version/series to choose expected record diffs. |
| sysrap/sstate.h | Repurposes m1group2 to hold both m1/m2 group velocities and adds accessors. |
| sysrap/sctx.h | Adds carried timing state (current_group_velocity, current_material_index). |
| src/g4app.h | Enables UseGivenVelocity(true) for optical photons during tracking. |
| qudarap/QState.hh | Initializes debug/mock sstate with both material group velocities. |
| qudarap/qsim.h | Uses carried group velocity for time advance; updates it on boundary outcomes. |
| qudarap/QSim.cu | Initializes new sctx fields in debug kernels. |
| qudarap/qbnd.h | Populates material2 group velocity into sstate during boundary lookup. |
💡 Add Copilot custom instructions for smarter, more guided reviews. Learn how to get started.
There was a problem hiding this comment.
💡 Codex Review
Here are some automated review suggestions for this pull request.
Reviewed commit: 7b14deab52
ℹ️ About Codex in GitHub
Your team has set up Codex to review pull requests in this repo. Reviews are triggered when you
- Open a pull request for review
- Mark a draft as ready
- Comment "@codex review".
If Codex has suggestions, it will comment; otherwise it will react with 👍.
Codex can also answer questions or update the PR. Try commenting "@codex address that feedback".
Co-authored-by: Copilot <175728472+Copilot@users.noreply.github.com>
Cpp-Linter Report
|
There was a problem hiding this comment.
Cpp-linter Review
Used clang-format v20.1.2
Only 6 out of 8 clang-format concerns fit within this pull request's diff.
Click here for the full clang-format patch
diff --git a/qudarap/QState.hh b/qudarap/QState.hh
index 8484bf6..fc08f27 100644
--- a/qudarap/QState.hh
+++ b/qudarap/QState.hh
@@ -53 +53 @@ inline sstate QState::Make()
- float m1_group_velocity = 300.f ;
+ float m1_group_velocity = 300.f ;
@@ -68,2 +68,2 @@ inline sstate QState::Make()
- s.material2 = make_float4( m2_refractive_index, m2_absorption_length, m2_scattering_length, m2_reemission_prob );
- s.m1group2 = make_float4(m1_group_velocity, m2_group_velocity, 0.f, 0.f);
+ s.material2 = make_float4( m2_refractive_index, m2_absorption_length, m2_scattering_length, m2_reemission_prob );
+ s.m1group2 = make_float4(m1_group_velocity, m2_group_velocity, 0.f, 0.f);
diff --git a/qudarap/qbnd.h b/qudarap/qbnd.h
index 57ff4a0..3b05ab7 100644
--- a/qudarap/qbnd.h
+++ b/qudarap/qbnd.h
@@ -127 +126,0 @@ inline QBND_METHOD float4 qbnd::boundary_lookup( float nm, unsigned line, unsign
-
@@ -196 +195,2 @@ inline QBND_METHOD void qbnd::fill_state(sstate& s, unsigned boundary, float wav
- s.m1group2 = boundary_lookup( wavelength, m1_line, 1); // x: material1 group_velocity, y: material2 group_velocity, z/w unused
+ s.m1group2 =
+ boundary_lookup(wavelength, m1_line, 1); // x: material1 group_velocity, y: material2 group_velocity, z/w unused
diff --git a/qudarap/qsim.h b/qudarap/qsim.h
index c5006f9..11f96e5 100644
--- a/qudarap/qsim.h
+++ b/qudarap/qsim.h
@@ -726 +726 @@ inline QSIM_METHOD int qsim::propagate_to_boundary(unsigned& flag, RNG& rng, sct
- const float& group_velocity = ctx.current_group_velocity;
+ const float &group_velocity = ctx.current_group_velocity;
diff --git a/sysrap/sstate.h b/sysrap/sstate.h
index dc6da0e..7b6acc5 100644
--- a/sysrap/sstate.h
+++ b/sysrap/sstate.h
@@ -28 +28 @@ struct sstate
- float4 m1group2 ; // material1_group_velocity/material2_group_velocity/spare2/spare3
+ float4 m1group2; // material1_group_velocity/material2_group_velocity/spare2/spare3
@@ -45,4 +45,12 @@ struct sstate
- SSTATE_METHOD float material1_group_velocity() const { return m1group2.x; }
- SSTATE_METHOD float material2_group_velocity() const { return m1group2.y; }
- SSTATE_METHOD void set_material2_group_velocity(float gv) { m1group2.y = gv; }
-
+ SSTATE_METHOD float material1_group_velocity() const
+ {
+ return m1group2.x;
+ }
+ SSTATE_METHOD float material2_group_velocity() const
+ {
+ return m1group2.y;
+ }
+ SSTATE_METHOD void set_material2_group_velocity(float gv)
+ {
+ m1group2.y = gv;
+ }
@@ -71,21 +79,10 @@ inline std::ostream& operator<<(std::ostream& os, const sstate& s )
- os << "sstate"
- << std::endl
- << " material1 " << s.material1
- << " (refractive_index/absorption_length/scattering_length/reemission_prob) "
- << std::endl
- << " m1group2 " << s.m1group2
- << " (material1_group_velocity/material2_group_velocity/spare2/spare3) "
- << std::endl
- << " material2 " << s.material2
- << " (refractive_index/absorption_length/scattering_length/reemission_prob) "
- << std::endl
- << " surface " << s.surface
- << " (detect/absorb/reflect_specular/reflect_diffuse) "
- << std::endl
- << " optical " << s.optical
- << " (x/y/z/w index/type/finish/value) "
- << std::endl
- << " index " << s.index
- << " (indices of m1/m2/surf/sensor) "
- << std::endl
- ;
+ os << "sstate" << std::endl
+ << " material1 " << s.material1 << " (refractive_index/absorption_length/scattering_length/reemission_prob) "
+ << std::endl
+ << " m1group2 " << s.m1group2 << " (material1_group_velocity/material2_group_velocity/spare2/spare3) "
+ << std::endl
+ << " material2 " << s.material2 << " (refractive_index/absorption_length/scattering_length/reemission_prob) "
+ << std::endl
+ << " surface " << s.surface << " (detect/absorb/reflect_specular/reflect_diffuse) " << std::endl
+ << " optical " << s.optical << " (x/y/z/w index/type/finish/value) " << std::endl
+ << " index " << s.index << " (indices of m1/m2/surf/sensor) " << std::endl;
Have any feedback or feature suggestions? Share it here.
Geant4 11.4 substantially refactored G4OpBoundaryProcess, including the
dielectric-dielectric boundary handling path. In our simg4ox and
GPURaytrace workflows, Geant4 11.4+ is the first release family where
timing-sensitive optical comparisons diverge from the Geant4 11.3
baselines. That makes post-boundary timing correctness explicit in a way
our previous implementation did not handle.
This branch fixes two gaps in our handling of post-boundary optical
timing:
The Opticks GPU path advanced the next segment using the pre-boundary
medium group velocity even after transmission into material2.
The lightweight Geant4-side tracking/recording path used by simg4ox
did not explicitly enable UseGivenVelocity(true) for optical photons,
so the track was not forced to follow the boundary-selected velocity
during post-boundary recording.
Fix the GPU side by carrying the active-medium timing state through qsim:
reflection/transmission
boundary orientation still matches the active medium
Fix the Geant4-side CPU recording path by enabling
UseGivenVelocity(true) for optical photons in g4app tracking action, so
the recorded post-boundary timing follows the velocity selected by
Geant4.
Because hit merging in this codebase is time-bucket based, correcting
post-boundary timing changes timing-sensitive observables and can change
merged hit counts. The follow-up test updates on this branch therefore
make expected residual diffs and hit-count baselines depend on the
Geant4 release family instead of treating the new values as regressions.