Geant4  10.02.p03
RE05StackingAction Class Reference

#include <RE05StackingAction.hh>

Inheritance diagram for RE05StackingAction:
Collaboration diagram for RE05StackingAction:

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)
 

Private Member Functions

G4bool InsideRoI (const G4Track *aTrack, G4double ang)
 
G4VHitsCollectionGetCollection (G4String colName)
 

Private Attributes

RE05TrackerHitsCollectiontrkHits
 
RE05MuonHitsCollectionmuonHits
 
RE05StackingActionMessengertheMessenger
 
G4int stage
 
G4int reqMuon
 
G4int reqIsoMuon
 
G4int reqIso
 
G4double angRoI
 

Additional Inherited Members

- Protected Attributes inherited from G4UserStackingAction
G4StackManagerstackManager
 

Detailed Description

Definition at line 45 of file RE05StackingAction.hh.

Constructor & Destructor Documentation

◆ RE05StackingAction()

RE05StackingAction::RE05StackingAction ( )

Definition at line 45 of file RE05StackingAction.cc.

46  : trkHits(0), muonHits(0), stage(0)
47 {
48  angRoI = 30.0*deg;
49  reqMuon = 2;
50  reqIso = 10;
52 }
RE05StackingActionMessenger * theMessenger
RE05MuonHitsCollection * muonHits
static const double deg
Definition: G4SIunits.hh:151
RE05TrackerHitsCollection * trkHits

◆ ~RE05StackingAction()

RE05StackingAction::~RE05StackingAction ( )
virtual

Definition at line 54 of file RE05StackingAction.cc.

55 { delete theMessenger; }
RE05StackingActionMessenger * theMessenger

Member Function Documentation

◆ ClassifyNewTrack()

G4ClassificationOfNewTrack RE05StackingAction::ClassifyNewTrack ( const G4Track *  aTrack)
virtual

Reimplemented from G4UserStackingAction.

Definition at line 58 of file RE05StackingAction.cc.

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 }
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:94
G4bool InsideRoI(const G4Track *aTrack, G4double ang)
static G4MuonMinus * MuonMinusDefinition()
Definition: G4MuonMinus.cc:95
Here is the call graph for this function:

◆ GetCollection()

G4VHitsCollection * RE05StackingAction::GetCollection ( G4String  colName)
private

Definition at line 120 of file RE05StackingAction.cc.

121 {
124  int colID = SDMan->GetCollectionID(colName);
125  if(colID>=0)
126  {
127  const G4Event* currentEvent = runMan->GetCurrentEvent();
128  G4HCofThisEvent* HCE = currentEvent->GetHCofThisEvent();
129  return HCE->GetHC(colID);
130  }
131  return 0;
132 }
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:135
const G4Event * GetCurrentEvent() const
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
G4HCofThisEvent * GetHCofThisEvent() const
Definition: G4Event.hh:185
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetNIsolation()

G4int RE05StackingAction::GetNIsolation ( ) const
inline

Definition at line 76 of file RE05StackingAction.hh.

76 { return reqIso; }
Here is the caller graph for this function:

◆ GetNRequestIsoMuon()

G4int RE05StackingAction::GetNRequestIsoMuon ( ) const
inline

Definition at line 74 of file RE05StackingAction.hh.

74 { return reqIsoMuon; }
Here is the caller graph for this function:

◆ GetNRequestMuon()

G4int RE05StackingAction::GetNRequestMuon ( ) const
inline

Definition at line 72 of file RE05StackingAction.hh.

72 { return reqMuon; }
Here is the caller graph for this function:

◆ GetRoIAngle()

G4double RE05StackingAction::GetRoIAngle ( ) const
inline

Definition at line 78 of file RE05StackingAction.hh.

78 { return angRoI; }
Here is the caller graph for this function:

◆ InsideRoI()

G4bool RE05StackingAction::InsideRoI ( const G4Track *  aTrack,
G4double  ang 
)
private

Definition at line 99 of file RE05StackingAction.cc.

100 {
101  if(!muonHits)
102  { muonHits = (RE05MuonHitsCollection*)GetCollection("muonCollection"); }
103  if(!muonHits)
104  { G4cerr << "muonCollection NOT FOUND" << G4endl;
105  return true; }
106 
107  G4int nhits = muonHits->entries();
108 
109  const G4ThreeVector trPos = aTrack->GetPosition();
110  for(G4int i=0;i<nhits;i++)
111  {
112  G4ThreeVector muHitPos = (*muonHits)[i]->GetPos();
113  G4double angl = muHitPos.angle(trPos);
114  if(angl<ang) { return true; }
115  }
116 
117  return false;
118 }
G4VHitsCollection * GetCollection(G4String colName)
RE05MuonHitsCollection * muonHits
int G4int
Definition: G4Types.hh:78
double angle(const Hep3Vector &) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4GLOB_DLL std::ostream G4cerr
Here is the call graph for this function:
Here is the caller graph for this function:

◆ NewStage()

void RE05StackingAction::NewStage ( )
virtual

Reimplemented from G4UserStackingAction.

Definition at line 134 of file RE05StackingAction.cc.

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();
150  if(nhits<reqMuon)
151  {
152  stackManager->clear();
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  }
188  if(isoMuon<reqIsoMuon)
189  {
190  stackManager->clear();
192  return;
193  }
195  return;
196  }
197 
198  else
199  {
200  // Other stage change : just re-classify
202  }
203 }
G4VHitsCollection * GetCollection(G4String colName)
RE05MuonHitsCollection * muonHits
int G4int
Definition: G4Types.hh:78
double angle(const Hep3Vector &) const
RE05TrackerHitsCollection * trkHits
#define G4endl
Definition: G4ios.hh:61
G4StackManager * stackManager
G4GLOB_DLL std::ostream G4cerr
Here is the call graph for this function:

◆ PrepareNewEvent()

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 }
RE05MuonHitsCollection * muonHits
RE05TrackerHitsCollection * trkHits

◆ SetNIsolation()

void RE05StackingAction::SetNIsolation ( G4int  val)
inline

Definition at line 75 of file RE05StackingAction.hh.

75 { reqIso = val; }
Here is the caller graph for this function:

◆ SetNRequestIsoMuon()

void RE05StackingAction::SetNRequestIsoMuon ( G4int  val)
inline

Definition at line 73 of file RE05StackingAction.hh.

73 { reqIsoMuon = val; }
Here is the caller graph for this function:

◆ SetNRequestMuon()

void RE05StackingAction::SetNRequestMuon ( G4int  val)
inline

Definition at line 71 of file RE05StackingAction.hh.

71 { reqMuon = val; }
Here is the caller graph for this function:

◆ SetRoIAngle()

void RE05StackingAction::SetRoIAngle ( G4double  val)
inline

Definition at line 77 of file RE05StackingAction.hh.

77 { angRoI = val; }
Here is the caller graph for this function:

Member Data Documentation

◆ angRoI

G4double RE05StackingAction::angRoI
private

Definition at line 68 of file RE05StackingAction.hh.

◆ muonHits

RE05MuonHitsCollection* RE05StackingAction::muonHits
private

Definition at line 61 of file RE05StackingAction.hh.

◆ reqIso

G4int RE05StackingAction::reqIso
private

Definition at line 67 of file RE05StackingAction.hh.

◆ reqIsoMuon

G4int RE05StackingAction::reqIsoMuon
private

Definition at line 66 of file RE05StackingAction.hh.

◆ reqMuon

G4int RE05StackingAction::reqMuon
private

Definition at line 65 of file RE05StackingAction.hh.

◆ stage

G4int RE05StackingAction::stage
private

Definition at line 64 of file RE05StackingAction.hh.

◆ theMessenger

RE05StackingActionMessenger* RE05StackingAction::theMessenger
private

Definition at line 62 of file RE05StackingAction.hh.

◆ trkHits

RE05TrackerHitsCollection* RE05StackingAction::trkHits
private

Definition at line 60 of file RE05StackingAction.hh.


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