Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
RE05StackingAction Class Reference

#include <RE05StackingAction.hh>

Inheritance diagram for RE05StackingAction:
G4UserStackingAction

Public Member Functions

 RE05StackingAction ()
 
virtual ~RE05StackingAction ()
 
virtual G4ClassificationOfNewTrack ClassifyNewTrack (const G4Track *aTrack)
 
virtual void NewStage ()
 
virtual void PrepareNewEvent ()
 
void SetNRequestMuon (G4int val)
 
G4int GetNRequestMuon () const
 
void SetNRequestIsoMuon (G4int val)
 
G4int GetNRequestIsoMuon () const
 
void SetNIsolation (G4int val)
 
G4int GetNIsolation () const
 
void SetRoIAngle (G4double val)
 
G4double GetRoIAngle () const
 
- Public Member Functions inherited from G4UserStackingAction
 G4UserStackingAction ()
 
virtual ~G4UserStackingAction ()
 
void SetStackManager (G4StackManager *value)
 

Additional Inherited Members

- Protected Attributes inherited from G4UserStackingAction
G4StackManagerstackManager
 

Detailed Description

Definition at line 45 of file RE05StackingAction.hh.

Constructor & Destructor Documentation

RE05StackingAction::RE05StackingAction ( )

Definition at line 45 of file RE05StackingAction.cc.

References python.hepunit::deg.

46  : trkHits(0), muonHits(0), stage(0)
47 {
48  angRoI = 30.0*deg;
49  reqMuon = 2;
50  reqIso = 10;
51  theMessenger = new RE05StackingActionMessenger(this);
52 }
RE05StackingAction::~RE05StackingAction ( )
virtual

Definition at line 54 of file RE05StackingAction.cc.

55 { delete theMessenger; }

Member Function Documentation

G4ClassificationOfNewTrack RE05StackingAction::ClassifyNewTrack ( const G4Track aTrack)
virtual

Reimplemented from G4UserStackingAction.

Definition at line 58 of file RE05StackingAction.cc.

References fKill, fSuspend, fUrgent, fWaiting, G4Track::GetDefinition(), G4Track::GetParentID(), G4ParticleDefinition::GetPDGCharge(), G4Track::GetTrackStatus(), G4MuonMinus::MuonMinusDefinition(), and G4MuonPlus::MuonPlusDefinition().

59 {
60  G4ClassificationOfNewTrack classification = fWaiting;
61  switch(stage)
62  {
63  case 0: // Stage 0 : Primary muons only
64  if(aTrack->GetParentID()==0)
65  {
66  G4ParticleDefinition * particleType = aTrack->GetDefinition();
67  if((particleType==G4MuonPlus::MuonPlusDefinition())
68  ||(particleType==G4MuonMinus::MuonMinusDefinition()))
69  { classification = fUrgent; }
70  }
71  break;
72 
73  case 1: // Stage 1 : Charged primaries only
74  // Suspended tracks will be sent to the waiting stack
75  if(aTrack->GetParentID()!=0) { break; }
76  if(aTrack->GetTrackStatus()==fSuspend) { break; }
77  if(aTrack->GetDefinition()->GetPDGCharge()==0.) { break; }
78  classification = fUrgent;
79  break;
80 
81  default: // Stage 2 : Accept all primaries
82  // Accept all secondaries in RoI
83  // Kill secondaries outside RoI
84  if(aTrack->GetParentID()==0)
85  {
86  classification = fUrgent;
87  break;
88  }
89  if((angRoI<0.)||InsideRoI(aTrack,angRoI))
90  {
91  classification = fUrgent;
92  break;
93  }
94  classification = fKill;
95  }
96  return classification;
97 }
G4ParticleDefinition * GetDefinition() const
G4int GetParentID() const
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:94
G4TrackStatus GetTrackStatus() const
static G4MuonMinus * MuonMinusDefinition()
Definition: G4MuonMinus.cc:95
G4double GetPDGCharge() const
G4int RE05StackingAction::GetNIsolation ( ) const
inline

Definition at line 76 of file RE05StackingAction.hh.

Referenced by RE05StackingActionMessenger::GetCurrentValue().

76 { return reqIso; }
G4int RE05StackingAction::GetNRequestIsoMuon ( ) const
inline

Definition at line 74 of file RE05StackingAction.hh.

Referenced by RE05StackingActionMessenger::GetCurrentValue().

74 { return reqIsoMuon; }
G4int RE05StackingAction::GetNRequestMuon ( ) const
inline

Definition at line 72 of file RE05StackingAction.hh.

Referenced by RE05StackingActionMessenger::GetCurrentValue().

72 { return reqMuon; }
G4double RE05StackingAction::GetRoIAngle ( ) const
inline

Definition at line 78 of file RE05StackingAction.hh.

Referenced by RE05StackingActionMessenger::GetCurrentValue().

78 { return angRoI; }
void RE05StackingAction::NewStage ( )
virtual

Reimplemented from G4UserStackingAction.

Definition at line 134 of file RE05StackingAction.cc.

References CLHEP::Hep3Vector::angle(), G4StackManager::clear(), G4THitsCollection< T >::entries(), G4cerr, G4endl, G4StackManager::ReClassify(), and G4UserStackingAction::stackManager.

135 {
136  stage++;
137  G4int nhits;
138  if(stage==1)
139  {
140  // Stage 0->1 : check if at least "reqMuon" hits on muon chamber
141  // otherwise abort current event
142  if(!muonHits)
143  { muonHits = (RE05MuonHitsCollection*)GetCollection("muonCollection"); }
144  if(!muonHits)
145  { G4cerr << "muonCollection NOT FOUND" << G4endl;
146  return; }
147  nhits = muonHits->entries();
148 //// G4cout << "Stage 0->1 : " << nhits << " hits found in the muon chamber."
149 //// << G4endl;
150  if(nhits<reqMuon)
151  {
152  stackManager->clear();
153 //// G4cout << "++++++++ event aborted" << G4endl;
154  return;
155  }
157  return;
158  }
159 
160  else if(stage==2)
161  {
162  // Stage 1->2 : check the isolation of muon tracks
163  // at least "reqIsoMuon" isolated muons
164  // otherwise abort current event.
165  // Isolation requires "reqIso" or less hits
166  // (including own hits) in the RoI region
167  // in the tracker layers.
168  nhits = muonHits->entries();
169  if(!trkHits)
170  { trkHits = (RE05TrackerHitsCollection*)GetCollection("trackerCollection"); }
171  if(!trkHits)
172  { G4cerr << "trackerCollection NOT FOUND" << G4endl;
173  return; }
174  G4int nTrkhits = trkHits->entries();
175  G4int isoMuon = 0;
176  for(G4int j=0;j<nhits;j++)
177  {
178  G4ThreeVector hitPos = (*muonHits)[j]->GetPos();
179  G4int nhitIn = 0;
180  for(G4int jj=0;(jj<nTrkhits)&&(nhitIn<=reqIso);jj++)
181  {
182  G4ThreeVector trkhitPos = (*trkHits)[jj]->GetPos();
183  if(trkhitPos.angle(hitPos)<angRoI) nhitIn++;
184  }
185  if(nhitIn<=reqIso) isoMuon++;
186  }
187 //// G4cout << "Stage 1->2 : " << isoMuon << " isolated muon found." << G4endl;
188  if(isoMuon<reqIsoMuon)
189  {
190  stackManager->clear();
191 //// G4cout << "++++++++ event aborted" << G4endl;
192  return;
193  }
195  return;
196  }
197 
198  else
199  {
200  // Other stage change : just re-classify
202  }
203 }
double angle(const Hep3Vector &) const
int G4int
Definition: G4Types.hh:78
#define G4endl
Definition: G4ios.hh:61
G4StackManager * stackManager
G4GLOB_DLL std::ostream G4cerr
void RE05StackingAction::PrepareNewEvent ( )
virtual

Reimplemented from G4UserStackingAction.

Definition at line 205 of file RE05StackingAction.cc.

206 {
207  stage = 0;
208  trkHits = 0;
209  muonHits = 0;
210 }
void RE05StackingAction::SetNIsolation ( G4int  val)
inline

Definition at line 75 of file RE05StackingAction.hh.

Referenced by RE05StackingActionMessenger::SetNewValue().

75 { reqIso = val; }
void RE05StackingAction::SetNRequestIsoMuon ( G4int  val)
inline

Definition at line 73 of file RE05StackingAction.hh.

Referenced by RE05StackingActionMessenger::SetNewValue().

73 { reqIsoMuon = val; }
void RE05StackingAction::SetNRequestMuon ( G4int  val)
inline

Definition at line 71 of file RE05StackingAction.hh.

Referenced by RE05StackingActionMessenger::SetNewValue().

71 { reqMuon = val; }
void RE05StackingAction::SetRoIAngle ( G4double  val)
inline

Definition at line 77 of file RE05StackingAction.hh.

Referenced by RE05StackingActionMessenger::SetNewValue().

77 { angRoI = val; }

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