Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RE05StackingAction.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: RE05StackingAction.cc 98775 2016-08-09 14:30:39Z gcosmo $
27 //
30 //
31 
32 #include "RE05StackingAction.hh"
33 #include "G4SDManager.hh"
34 #include "G4RunManager.hh"
35 #include "G4Event.hh"
36 #include "G4HCofThisEvent.hh"
37 #include "G4Track.hh"
38 #include "G4TrackStatus.hh"
39 #include "G4ParticleDefinition.hh"
40 #include "G4ParticleTypes.hh"
42 #include "G4SystemOfUnits.hh"
43 #include "G4ios.hh"
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
46 
49  fTrkHits(0), fMuonHits(0), fMessenger(0),
50  fStage(0),
51  fReqMuon(2),
52  fReqIso(10),
53  fAngRoI(30.0*deg)
54 {
55  fMessenger = new RE05StackingActionMessenger(this);
56 }
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59 
61 { delete fMessenger; }
62 
63 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64 
67 {
68  G4ClassificationOfNewTrack classification = fWaiting;
69  switch(fStage)
70  {
71  case 0: // Stage 0 : Primary muons only
72  if(aTrack->GetParentID()==0)
73  {
74  G4ParticleDefinition * particleType = aTrack->GetDefinition();
75  if((particleType==G4MuonPlus::MuonPlusDefinition())
76  ||(particleType==G4MuonMinus::MuonMinusDefinition()))
77  { classification = fUrgent; }
78  }
79  break;
80 
81  case 1: // Stage 1 : Charged primaries only
82  // Suspended tracks will be sent to the waiting stack
83  if(aTrack->GetParentID()!=0) { break; }
84  if(aTrack->GetTrackStatus()==fSuspend) { break; }
85  if(aTrack->GetDefinition()->GetPDGCharge()==0.) { break; }
86  classification = fUrgent;
87  break;
88 
89  default: // Stage 2 : Accept all primaries
90  // Accept all secondaries in RoI
91  // Kill secondaries outside RoI
92  if(aTrack->GetParentID()==0)
93  {
94  classification = fUrgent;
95  break;
96  }
97  if((fAngRoI<0.)||InsideRoI(aTrack,fAngRoI))
98  {
99  classification = fUrgent;
100  break;
101  }
102  classification = fKill;
103  }
104  return classification;
105 }
106 
107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108 
109 G4bool RE05StackingAction::InsideRoI(const G4Track * aTrack,G4double ang)
110 {
111  if(!fMuonHits)
112  { fMuonHits = (RE05MuonHitsCollection*)GetCollection("muonCollection"); }
113  if(!fMuonHits)
114  { G4cerr << "muonCollection NOT FOUND" << G4endl;
115  return true; }
116 
117  G4int nhits = fMuonHits->entries();
118 
119  const G4ThreeVector trPos = aTrack->GetPosition();
120  for(G4int i=0;i<nhits;i++)
121  {
122  G4ThreeVector muHitPos = (*fMuonHits)[i]->GetPos();
123  G4double angl = muHitPos.angle(trPos);
124  if(angl<ang) { return true; }
125  }
126 
127  return false;
128 }
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
131 
132 G4VHitsCollection* RE05StackingAction::GetCollection(G4String colName)
133 {
136  int colID = SDMan->GetCollectionID(colName);
137  if(colID>=0)
138  {
139  const G4Event* currentEvent = runMan->GetCurrentEvent();
140  G4HCofThisEvent* HCE = currentEvent->GetHCofThisEvent();
141  return HCE->GetHC(colID);
142  }
143  return 0;
144 }
145 
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
147 
149 {
150  fStage++;
151  G4int nhits;
152  if(fStage==1)
153  {
154  // Stage 0->1 : check if at least "fReqMuon" hits on muon chamber
155  // otherwise abort current event
156  if(!fMuonHits)
157  { fMuonHits = (RE05MuonHitsCollection*)GetCollection("muonCollection"); }
158  if(!fMuonHits)
159  { G4cerr << "muonCollection NOT FOUND" << G4endl;
160  return; }
161  nhits = fMuonHits->entries();
164  if(nhits<fReqMuon)
165  {
166  stackManager->clear();
168  return;
169  }
171  return;
172  }
173 
174  else if(fStage==2)
175  {
176  // Stage 1->2 : check the isolation of muon tracks
177  // at least "fReqIsoMuon" isolated muons
178  // otherwise abort current event.
179  // Isolation requires "fReqIso" or less hits
180  // (including own hits) in the RoI region
181  // in the tracker layers.
182  nhits = fMuonHits->entries();
183  if(!fTrkHits)
184  { fTrkHits = (RE05TrackerHitsCollection*)GetCollection("trackerCollection"); }
185  if(!fTrkHits)
186  { G4cerr << "trackerCollection NOT FOUND" << G4endl;
187  return; }
188  G4int nTrkhits = fTrkHits->entries();
189  G4int isoMuon = 0;
190  for(G4int j=0;j<nhits;j++)
191  {
192  G4ThreeVector hitPos = (*fMuonHits)[j]->GetPos();
193  G4int nhitIn = 0;
194  for(G4int jj=0;(jj<nTrkhits)&&(nhitIn<=fReqIso);jj++)
195  {
196  G4ThreeVector trkhitPos = (*fTrkHits)[jj]->GetPos();
197  if(trkhitPos.angle(hitPos)<fAngRoI) nhitIn++;
198  }
199  if(nhitIn<=fReqIso) isoMuon++;
200  }
202  if(isoMuon<fReqIsoMuon)
203  {
204  stackManager->clear();
206  return;
207  }
209  return;
210  }
211 
212  else
213  {
214  // Other fStage change : just re-classify
216  }
217 }
218 
219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
220 
222 {
223  fStage = 0;
224  fTrkHits = 0;
225  fMuonHits = 0;
226 }
227 
228 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4ParticleDefinition * GetDefinition() const
G4int GetParentID() const
static G4MuonPlus * MuonPlusDefinition()
Definition: G4MuonPlus.cc:94
G4VHitsCollection * GetHC(G4int i)
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:135
double angle(const Hep3Vector &) const
Definition of the RE05StackingActionMessenger class.
const G4ThreeVector & GetPosition() const
G4TrackStatus GetTrackStatus() const
G4int entries() const
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
Definition of the RE05StackingAction class.
virtual G4ClassificationOfNewTrack ClassifyNewTrack(const G4Track *aTrack)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
static G4MuonMinus * MuonMinusDefinition()
Definition: G4MuonMinus.cc:95
virtual void PrepareNewEvent()
#define G4endl
Definition: G4ios.hh:61
G4HCofThisEvent * GetHCofThisEvent() const
Definition: G4Event.hh:185
G4StackManager * stackManager
double G4double
Definition: G4Types.hh:76
static constexpr double deg
Definition: G4SIunits.hh:152
const G4Event * GetCurrentEvent() const
G4double GetPDGCharge() const
G4GLOB_DLL std::ostream G4cerr