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

Event action. More...

#include <B5EventAction.hh>

Inheritance diagram for B5EventAction:
G4UserEventAction

Public Member Functions

 B5EventAction ()
 
virtual ~B5EventAction ()
 
virtual void BeginOfEventAction (const G4Event *)
 
virtual void EndOfEventAction (const G4Event *)
 
- Public Member Functions inherited from G4UserEventAction
 G4UserEventAction ()
 
virtual ~G4UserEventAction ()
 
void SetEventManager (G4EventManager *value)
 

Additional Inherited Members

- Protected Attributes inherited from G4UserEventAction
G4EventManagerfpEventManager
 

Detailed Description

Event action.

Definition at line 40 of file B5EventAction.hh.

Constructor & Destructor Documentation

B5EventAction::B5EventAction ( )

Definition at line 49 of file B5EventAction.cc.

References G4RunManager::GetRunManager(), and G4RunManager::SetPrintProgress().

51  fHHC1ID(-1),
52  fHHC2ID(-1),
53  fDHC1ID(-1),
54  fDHC2ID(-1),
55  fECHCID(-1),
56  fHCHCID(-1)
57 {
58  // set printing per each event
60 }
void SetPrintProgress(G4int i)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74
B5EventAction::~B5EventAction ( )
virtual

Definition at line 64 of file B5EventAction.cc.

65 {}

Member Function Documentation

void B5EventAction::BeginOfEventAction ( const G4Event )
virtual

Reimplemented from G4UserEventAction.

Definition at line 69 of file B5EventAction.cc.

References G4SDManager::GetCollectionID(), and G4SDManager::GetSDMpointer().

70 {
71  if (fHHC1ID==-1) {
73  fHHC1ID = sdManager->GetCollectionID("hodoscope1/hodoscopeColl");
74  fHHC2ID = sdManager->GetCollectionID("hodoscope2/hodoscopeColl");
75  fDHC1ID = sdManager->GetCollectionID("chamber1/driftChamberColl");
76  fDHC2ID = sdManager->GetCollectionID("chamber2/driftChamberColl");
77  fECHCID = sdManager->GetCollectionID("EMcalorimeter/EMcalorimeterColl");
78  fHCHCID = sdManager->GetCollectionID("HadCalorimeter/HadCalorimeterColl");
79  }
80 }
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:131
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
void B5EventAction::EndOfEventAction ( const G4Event event)
virtual

Reimplemented from G4UserEventAction.

Definition at line 84 of file B5EventAction.cc.

References G4VAnalysisManager::AddNtupleRow(), G4THitsCollection< T >::entries(), G4VAnalysisManager::FillH1(), G4VAnalysisManager::FillH2(), G4VAnalysisManager::FillNtupleDColumn(), G4VAnalysisManager::FillNtupleIColumn(), G4cout, G4endl, G4Exception(), B5EmCalorimeterHit::GetEdep(), B5HadCalorimeterHit::GetEdep(), G4Event::GetEventID(), G4PrimaryParticle::GetG4code(), G4HCofThisEvent::GetHC(), B5DriftChamberHit::GetLayerID(), B5DriftChamberHit::GetLocalPos(), G4PrimaryParticle::GetMomentum(), G4ParticleDefinition::GetParticleName(), G4RunManager::GetPrintProgress(), G4RunManager::GetRunManager(), JustWarning, python.hepunit::MeV, B5HodoscopeHit::Print(), B5DriftChamberHit::Print(), CLHEP::Hep3Vector::x(), and CLHEP::Hep3Vector::y().

85 {
86  G4HCofThisEvent* hce = event->GetHCofThisEvent();
87  if (!hce)
88  {
90  msg << "No hits collection of this event found.\n";
91  G4Exception("B5EventAction::EndOfEventAction()",
92  "B5Code001", JustWarning, msg);
93  return;
94  }
95 
96 
97  // Get hits collections
99  = static_cast<B5HodoscopeHitsCollection*>(hce->GetHC(fHHC1ID));
100 
102  = static_cast<B5HodoscopeHitsCollection*>(hce->GetHC(fHHC2ID));
103 
105  = static_cast<B5DriftChamberHitsCollection*>(hce->GetHC(fDHC1ID));
106 
108  = static_cast<B5DriftChamberHitsCollection*>(hce->GetHC(fDHC2ID));
109 
111  = static_cast<B5EmCalorimeterHitsCollection*>(hce->GetHC(fECHCID));
112 
114  = static_cast<B5HadCalorimeterHitsCollection*>(hce->GetHC(fHCHCID));
115 
116  if ( (!hHC1) || (!hHC2) || (!dHC1) || (!dHC2) || (!ecHC) || (!hcHC) )
117  {
119  msg << "Some of hits collections of this event not found.\n";
120  G4Exception("B5EventAction::EndOfEventAction()",
121  "B5Code001", JustWarning, msg);
122  return;
123  }
124 
125  //
126  // Fill histograms & ntuple
127  //
128 
129  // Get analysis manager
130  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
131 
132  // Fill histograms
133 
134  G4int n_hit = dHC1->entries();
135  analysisManager->FillH1(0, n_hit);
136 
137  for (G4int i=0;i<n_hit;i++)
138  {
139  B5DriftChamberHit* hit = (*dHC1)[i];
140  G4ThreeVector localPos = hit->GetLocalPos();
141  analysisManager->FillH2(0, localPos.x(), localPos.y());
142  }
143 
144  n_hit = dHC2->entries();
145  analysisManager->FillH1(1, n_hit);
146 
147  for (G4int i=0;i<n_hit;i++)
148  {
149  B5DriftChamberHit* hit = (*dHC2)[i];
150  G4ThreeVector localPos = hit->GetLocalPos();
151  analysisManager->FillH2(1, localPos.x(), localPos.y());
152  }
153 
154 
155  // Fill ntuple
156 
157  // Dc1Hits
158  analysisManager->FillNtupleIColumn(0, dHC1->entries());
159  // Dc2Hits
160  analysisManager->FillNtupleIColumn(1, dHC1->entries());
161 
162  // ECEnergy
163  G4int totalEmHit = 0;
164  G4double totalEmE = 0.;
165  for (G4int i=0;i<80;i++)
166  {
167  B5EmCalorimeterHit* hit = (*ecHC)[i];
168  G4double eDep = hit->GetEdep();
169  if (eDep>0.)
170  {
171  totalEmHit++;
172  totalEmE += eDep;
173  }
174  }
175  analysisManager->FillNtupleDColumn(2, totalEmE);
176 
177  // HCEnergy
178  G4int totalHadHit = 0;
179  G4double totalHadE = 0.;
180  for (G4int i=0;i<20;i++)
181  {
182  B5HadCalorimeterHit* hit = (*hcHC)[i];
183  G4double eDep = hit->GetEdep();
184  if (eDep>0.)
185  {
186  totalHadHit++;
187  totalHadE += eDep;
188  }
189  }
190  analysisManager->FillNtupleDColumn(3, totalHadE);
191 
192  // Time 1
193  for (G4int i=0;i<hHC1->entries();i++)
194  {
195  analysisManager->FillNtupleDColumn(4,(*hHC1)[i]->GetTime());
196  }
197 
198  // Time 2
199  for (G4int i=0;i<hHC2->entries();i++)
200  {
201  analysisManager->FillNtupleDColumn(5,(*hHC2)[i]->GetTime());
202  }
203 
204  analysisManager->AddNtupleRow();
205 
206  //
207  // Print diagnostics
208  //
209 
211  if ( printModulo==0 || event->GetEventID() % printModulo != 0) return;
212 
213  G4PrimaryParticle* primary = event->GetPrimaryVertex(0)->GetPrimary(0);
214  G4cout << G4endl
215  << ">>> Event " << event->GetEventID() << " >>> Simulation truth : "
216  << primary->GetG4code()->GetParticleName()
217  << " " << primary->GetMomentum() << G4endl;
218 
219  // Hodoscope 1
220  n_hit = hHC1->entries();
221  G4cout << "Hodoscope 1 has " << n_hit << " hits." << G4endl;
222  for (G4int i=0;i<n_hit;i++)
223  {
224  B5HodoscopeHit* hit = (*hHC1)[i];
225  hit->Print();
226  }
227 
228  // Hodoscope 2
229  n_hit = hHC2->entries();
230  G4cout << "Hodoscope 2 has " << n_hit << " hits." << G4endl;
231  for (G4int i=0;i<n_hit;i++)
232  {
233  B5HodoscopeHit* hit = (*hHC2)[i];
234  hit->Print();
235  }
236 
237  // Drift chamber 1
238  n_hit = dHC1->entries();
239  G4cout << "Drift Chamber 1 has " << n_hit << " hits." << G4endl;
240  for (G4int i2=0;i2<5;i2++)
241  {
242  for (G4int i=0;i<n_hit;i++)
243  {
244  B5DriftChamberHit* hit = (*dHC1)[i];
245  if (hit->GetLayerID()==i2) hit->Print();
246  }
247  }
248 
249  // Drift chamber 2
250  n_hit = dHC2->entries();
251  G4cout << "Drift Chamber 2 has " << n_hit << " hits." << G4endl;
252  for (G4int i2=0;i2<5;i2++)
253  {
254  for (G4int i=0;i<n_hit;i++)
255  {
256  B5DriftChamberHit* hit = (*dHC2)[i];
257  if (hit->GetLayerID()==i2) hit->Print();
258  }
259  }
260 
261  // EM calorimeter
262  G4cout << "EM Calorimeter has " << totalEmHit << " hits. Total Edep is "
263  << totalEmE/MeV << " (MeV)" << G4endl;
264 
265  // Had calorimeter
266  G4cout << "Hadron Calorimeter has " << totalHadHit << " hits. Total Edep is "
267  << totalHadE/MeV << " (MeV)" << G4endl;
268 }
G4int GetPrintProgress()
G4ThreeVector GetMomentum() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
double x() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4ParticleDefinition * GetG4code() const
G4int GetEventID() const
Definition: G4Event.hh:140
G4ThreeVector GetLocalPos() const
G4bool FillNtupleIColumn(G4int id, G4int value)
G4GLOB_DLL std::ostream G4cout
G4double GetEdep() const
G4bool FillNtupleDColumn(G4int id, G4double value)
virtual void Print()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:74
G4bool FillH2(G4int id, G4double xvalue, G4double yvalue, G4double weight=1.0)
G4bool FillH1(G4int id, G4double value, G4double weight=1.0)
double y() const
G4double GetEdep() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4int GetLayerID() const

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