ATLAS Offline Software
Loading...
Searching...
No Matches
G4UA::AthenaStackingAction Class Reference

Standard ATLAS stacking action functionality. More...

#include <AthenaStackingAction.h>

Inheritance diagram for G4UA::AthenaStackingAction:
Collaboration diagram for G4UA::AthenaStackingAction:

Classes

struct  Config
 Configuration option struct for AthenaStackingAction. More...

Public Member Functions

 AthenaStackingAction (const Config &config)
 Constructor with configuration.
virtual G4ClassificationOfNewTrack ClassifyNewTrack (const G4Track *track) override
 Classify a new track.

Protected Member Functions

bool isNeutrino (const G4Track *) const
 Identify track as a neutrino.
bool isGamma (const G4Track *) const
 Identify track as a photon.
bool isNeutron (const G4Track *) const
 Identify track as a neutron.
PrimaryParticleInformationgetPrimaryParticleInformation (const G4Track *track) const
 obtain the PrimaryParticleInformation from the current G4Track

Protected Attributes

Config m_config
 Configuration options.
double m_oneOverWeightNeutron
double m_oneOverWeightPhoton

Detailed Description

Standard ATLAS stacking action functionality.

Mostly taken from the old AthenaStackingAction implementation, but adapted for the new user action design for multi-threading.

Author
Steve Farrell Steve.nosp@m.n.Fa.nosp@m.rrell.nosp@m.@cer.nosp@m.n.ch

Definition at line 22 of file AthenaStackingAction.h.

Constructor & Destructor Documentation

◆ AthenaStackingAction()

G4UA::AthenaStackingAction::AthenaStackingAction ( const Config & config)

Constructor with configuration.

Definition at line 42 of file AthenaStackingAction.cxx.

42 :
43 m_config(config),
46 {
47 // calculate this division only once
48 if (m_config.applyNRR)
49 m_oneOverWeightNeutron = 1./m_config.russianRouletteNeutronWeight;
50
51 // calculate this division only once
52 if (m_config.applyPRR)
53 m_oneOverWeightPhoton = 1./m_config.russianRoulettePhotonWeight;
54 }
Config m_config
Configuration options.

Member Function Documentation

◆ ClassifyNewTrack()

G4ClassificationOfNewTrack G4UA::AthenaStackingAction::ClassifyNewTrack ( const G4Track * track)
overridevirtual

Classify a new track.

Result can be fUrgent, fWaiting, fPostpone, or fKill.

Pass ownership to track. The G4VUserTrackInformation* fpUserInformation member variable set by this method is mutable. G4Tracks are thread-local.

Pass ownership to track. The G4VUserTrackInformation* fpUserInformation member variable set by this method is mutable. G4Tracks are thread-local.

Reimplemented in G4UA::AthenaDebugStackingAction.

Definition at line 60 of file AthenaStackingAction.cxx.

61 {
62 // Kill neutrinos if enabled
63 if(m_config.killAllNeutrinos && isNeutrino(track)) {
64 return fKill;
65 }
66
67 // Kill super-low-E photons
68 const double safeCut = 0.00005;
69 double totalE = track->GetTotalEnergy();
70 if(isGamma(track) && totalE < safeCut) {
71 return fKill;
72 }
73
74 // TODO: Why is this here? Can I remove it?
75 G4Event* ev = G4EventManager::GetEventManager()->GetNonconstCurrentEvent();
76 AtlasG4EventUserInfo* atlasG4EvtUserInfo __attribute__ ((unused)) =
77 static_cast<AtlasG4EventUserInfo*> (ev->GetUserInformation());
78
79 // Neutron Russian Roulette
80 if (m_config.applyNRR && isNeutron(track) &&
81 CxxUtils::fpcompare::equal(track->GetWeight(), 1.0) && // do not re-Roulette particles
82 track->GetKineticEnergy() < m_config.russianRouletteNeutronThreshold) {
83 // shoot random number
84 if ( CLHEP::RandFlat::shoot() > m_oneOverWeightNeutron ) {
85 // Kill (w-1)/w neutrons
86 return fKill;
87 }
88 // TODO There may be another way to set the weights via
89 // another G4 interface avoiding the const_cast, but the
90 // changes are more major and will need more careful validation.
91 G4Track* mutableTrack ATLAS_THREAD_SAFE = const_cast<G4Track*> (track);
92 // Weight the rest 1/w neutrons with a weight of w
93 mutableTrack->SetWeight(m_config.russianRouletteNeutronWeight);
94 }
95
96 // Photon Russian Roulette
97 if (m_config.applyPRR && isGamma(track) && track->GetOriginTouchable() &&
98 track->GetOriginTouchable()->GetVolume()->GetName().substr(0, 3) == "LAr" && // only for photons created in LAr
99 CxxUtils::fpcompare::equal(track->GetWeight(), 1.0) && // do not re-Roulette particles
100 track->GetKineticEnergy() < m_config.russianRoulettePhotonThreshold) {
101 // shoot random number
102 if ( CLHEP::RandFlat::shoot() > m_oneOverWeightPhoton ) {
103 // Kill (w-1)/w photons
104 return fKill;
105 }
106 // TODO There may be another way to set the weights via
107 // another G4 interface avoiding the const_cast, but the
108 // changes are more major and will need more careful validation.
109 G4Track* mutableTrack ATLAS_THREAD_SAFE = const_cast<G4Track*> (track);
110 // Weight the rest 1/w neutrons with a weight of w
111 mutableTrack->SetWeight(m_config.russianRoulettePhotonWeight);
112 }
113
114 // Handle primary particles
115 if(track->GetParentID() == 0) { // Condition for Primaries
116 // Extract the PrimaryParticleInformation
117 PrimaryParticleInformation* primaryPartInfo = this->getPrimaryParticleInformation(track);
118 // Fill some information for this track
119 if(primaryPartInfo) {
120 if (!m_config.isISFJob) {
121 // don't do anything
122 auto part = primaryPartInfo->GetHepMCParticle();
123 if (part) {
124 // OK, we got back to HepMC
125 std::unique_ptr<TrackInformation> ti = std::make_unique<TrackInformation>(part);
126 ti->SetRegenerationNr(0);
127 ti->SetClassification(TrackInformation::Primary);
128 // regNr=0 and classify=Primary are default values anyway
132 track->SetUserInformation(ti.release());
133 }
134 // TODO What does this condition mean?
135 else if(primaryPartInfo->GetParticleUniqueID() >= 0 && primaryPartInfo->GetParticleBarcode() >= 0) {
136 // PrimaryParticleInformation should at least provide a barcode
137 std::unique_ptr<TrackBarcodeInfo> bi = std::make_unique<TrackBarcodeInfo>(primaryPartInfo->GetParticleUniqueID(), primaryPartInfo->GetParticleBarcode());
141 track->SetUserInformation(bi.release());
142 }
143 } // no ISFParticle attached
144 } // has PrimaryParticleInformation
145 }
146 // Secondary track; decide whether to save or kill
147 else if( isGamma(track) &&
148 m_config.photonEnergyCut > 0 &&
149 totalE < m_config.photonEnergyCut )
150 {
151 return fKill;
152 }
153 return fUrgent;
154 }
__attribute__((always_inline)) inline uint16_t TileCalibDrawerBase
void unused(Args &&...)
#define ATLAS_THREAD_SAFE
bool isNeutrino(const G4Track *) const
Identify track as a neutrino.
bool isNeutron(const G4Track *) const
Identify track as a neutron.
bool isGamma(const G4Track *) const
Identify track as a photon.
PrimaryParticleInformation * getPrimaryParticleInformation(const G4Track *track) const
obtain the PrimaryParticleInformation from the current G4Track
HepMC::ConstGenParticlePtr GetHepMCParticle() const
return a pointer to the GenParticle used to create the G4PrimaryParticle
int ev
Definition globals.cxx:25
bool equal(double a, double b)
Compare two FP numbers, working around x87 precision issues.
Definition fpcompare.h:114

◆ getPrimaryParticleInformation()

PrimaryParticleInformation * G4UA::AthenaStackingAction::getPrimaryParticleInformation ( const G4Track * track) const
protected

obtain the PrimaryParticleInformation from the current G4Track

Definition at line 156 of file AthenaStackingAction.cxx.

157 {
158 const G4DynamicParticle* dp = track->GetDynamicParticle();
159 if(dp) {
160 const G4PrimaryParticle* pp = nullptr;
161 pp = dp->GetPrimaryParticle();
162 if(pp) {
163 // Extract the PrimaryParticleInformation
164 return dynamic_cast<PrimaryParticleInformation*>
165 ( pp->GetUserInformation() );
166 }
167 }
168 return nullptr;
169 }

◆ isGamma()

bool G4UA::AthenaStackingAction::isGamma ( const G4Track * track) const
protected

Identify track as a photon.

Definition at line 186 of file AthenaStackingAction.cxx.

187 {
188 return track->GetParticleDefinition() == G4Gamma::Gamma();
189 }

◆ isNeutrino()

bool G4UA::AthenaStackingAction::isNeutrino ( const G4Track * track) const
protected

Identify track as a neutrino.

It might be useful to move this kind of functionality into some standalong helper function(s).

Definition at line 174 of file AthenaStackingAction.cxx.

175 {
176 auto particleDef = track->GetParticleDefinition();
177 return (particleDef == G4NeutrinoE::NeutrinoEDefinition() ||
178 particleDef == G4AntiNeutrinoE::AntiNeutrinoEDefinition() ||
179 particleDef == G4NeutrinoMu::NeutrinoMuDefinition() ||
180 particleDef == G4AntiNeutrinoMu::AntiNeutrinoMuDefinition() ||
181 particleDef == G4NeutrinoTau::NeutrinoTauDefinition() ||
182 particleDef == G4AntiNeutrinoTau::AntiNeutrinoTauDefinition());
183 }

◆ isNeutron()

bool G4UA::AthenaStackingAction::isNeutron ( const G4Track * track) const
protected

Identify track as a neutron.

Definition at line 192 of file AthenaStackingAction.cxx.

193 {
194 return track->GetParticleDefinition() == G4Neutron::Neutron();
195 }

Member Data Documentation

◆ m_config

Config G4UA::AthenaStackingAction::m_config
protected

Configuration options.

Definition at line 61 of file AthenaStackingAction.h.

◆ m_oneOverWeightNeutron

double G4UA::AthenaStackingAction::m_oneOverWeightNeutron
protected

Definition at line 78 of file AthenaStackingAction.h.

◆ m_oneOverWeightPhoton

double G4UA::AthenaStackingAction::m_oneOverWeightPhoton
protected

Definition at line 81 of file AthenaStackingAction.h.


The documentation for this class was generated from the following files: