| Summary: | G4VParticleChange::SetNumberOfSecondaries() Warning theListOfSecondaries is not empty | ||
|---|---|---|---|
| Product: | Geant4 | Reporter: | Javier Muñoz <javier.munoz.vidal> |
| Component: | processes/optical | Assignee: | Daren Sawkey <daren.sawkey> |
| Status: | RESOLVED FIXED | ||
| Severity: | minor | ||
| Priority: | P4 | ||
| Version: | 10.5 | ||
| Hardware: | Apple | ||
| OS: | Mac OS X | ||
|
Description
Javier Muñoz
2019-10-24 19:41:53 CEST
Thanks for reporting this. It is now fixed for the upcoming release. Does it have a significant impact on the results got with current version ?? It can change the results, when the absorbed photon has an energy close to the low energy limit of the emission spectrum. (This might not be realistic.) For each emitted photon, the code will try to sample the energy from the spectrum, such that the emitted photon energy is less than the absorbed photon energy. If this can't be done, the number of secondaries is reduced by 1. This unintentionally deleted any photons that were previously generated this step.
The diff is below (changes to printouts suppressed). If you like, you can try it and report back.
diff --git a/source/processes/optical/src/G4OpWLS.cc b/source/processes/optical/src/G4OpWLS.cc
index 62469cdd82..639adfcc68 100644
--- a/source/processes/optical/src/G4OpWLS.cc
+++ b/source/processes/optical/src/G4OpWLS.cc
@@ -100,12 +100,13 @@ G4OpWLS::~G4OpWLS()
G4VParticleChange*
G4OpWLS::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
{
+ std::vector<G4Track*> proposedSecondaries;
aParticleChange.Initialize(aTrack);
aParticleChange.ProposeTrackStatus(fStopAndKill);
@@ -144,8 +145,6 @@ G4OpWLS::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
}
- aParticleChange.SetNumberOfSecondaries(NumPhotons);
-
G4double primaryEnergy = aTrack.GetDynamicParticle()->GetKineticEnergy();
G4int materialIndex = aMaterial->GetIndex();
@@ -191,19 +190,26 @@ G4OpWLS::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
// If no such energy can be sampled, return one less secondary, or none
if (sampledEnergy > primaryEnergy) {
NumberOfPhotons--;
- aParticleChange.SetNumberOfSecondaries(NumberOfPhotons);
if (NumberOfPhotons == 0) {
// return unchanged particle and no secondaries
+ aParticleChange.SetNumberOfSecondaries(0);
return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
}
continue;
@@ -268,10 +274,14 @@ G4OpWLS::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
// aSecondaryTrack->SetTouchableHandle((G4VTouchable*)0);
aSecondaryTrack->SetParentID(aTrack.GetTrackID());
-
- aParticleChange.AddSecondary(aSecondaryTrack);
+
+ proposedSecondaries.push_back(aSecondaryTrack);
}
+ aParticleChange.SetNumberOfSecondaries(proposedSecondaries.size());
+ for (auto sec : proposedSecondaries) {
+ aParticleChange.AddSecondary(sec);
+ }
if (verboseLevel>0) {
G4cout << "\n Exiting from G4OpWLS::DoIt -- NumberOfSecondaries = "
<< aParticleChange.GetNumberOfSecondaries() << G4endl;
|