Loading...
Searching...
No Matches
EventAction.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//
27/// \file B5/src/EventAction.cc
28/// \brief Implementation of the B5::EventAction class
29
30#include "EventAction.hh"
31#include "HodoscopeHit.hh"
32#include "DriftChamberHit.hh"
33#include "EmCalorimeterHit.hh"
34#include "HadCalorimeterHit.hh"
35#include "Constants.hh"
36
37#include "G4Event.hh"
38#include "G4RunManager.hh"
39#include "G4EventManager.hh"
40#include "G4HCofThisEvent.hh"
41#include "G4VHitsCollection.hh"
42#include "G4SDManager.hh"
43#include "G4SystemOfUnits.hh"
44#include "G4ios.hh"
45#include "G4AnalysisManager.hh"
46
47using std::array;
48using std::vector;
49
50namespace {
51
52// Utility function which finds a hit collection with the given Id
53// and print warnings if not found
54G4VHitsCollection* GetHC(const G4Event* event, G4int collId) {
55 auto hce = event->GetHCofThisEvent();
56 if (!hce) {
57 G4ExceptionDescription msg;
58 msg << "No hits collection of this event found." << G4endl;
59 G4Exception("EventAction::EndOfEventAction()",
60 "Code001", JustWarning, msg);
61 return nullptr;
62 }
63
64 auto hc = hce->GetHC(collId);
65 if ( ! hc) {
66 G4ExceptionDescription msg;
67 msg << "Hits collection " << collId << " of this event not found." << G4endl;
68 G4Exception("EventAction::EndOfEventAction()",
69 "Code001", JustWarning, msg);
70 }
71 return hc;
72}
73
74}
75
76namespace B5
77{
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80
82{
83 // set printing per each event
84 G4RunManager::GetRunManager()->SetPrintProgress(1);
85}
86
87//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
88
90{
91 // Find hit collections and histogram Ids by names (just once)
92 // and save them in the data members of this class
93
94 if (fHodHCID[0] == -1) {
95 auto sdManager = G4SDManager::GetSDMpointer();
96 auto analysisManager = G4AnalysisManager::Instance();
97
98 // hits collections names
99 array<G4String, kDim> hHCName
100 = {{ "hodoscope1/hodoscopeColl", "hodoscope2/hodoscopeColl" }};
101 array<G4String, kDim> dHCName
102 = {{ "chamber1/driftChamberColl", "chamber2/driftChamberColl" }};
103 array<G4String, kDim> cHCName
104 = {{ "EMcalorimeter/EMcalorimeterColl", "HadCalorimeter/HadCalorimeterColl" }};
105
106 // histograms names
107 array<array<G4String, kDim>, kDim> histoName
108 = {{ {{ "Chamber1", "Chamber2" }}, {{ "Chamber1 XY", "Chamber2 XY" }} }};
109
110 for (G4int iDet = 0; iDet < kDim; ++iDet) {
111 // hit collections IDs
112 fHodHCID[iDet] = sdManager->GetCollectionID(hHCName[iDet]);
113 fDriftHCID[iDet] = sdManager->GetCollectionID(dHCName[iDet]);
114 fCalHCID[iDet] = sdManager->GetCollectionID(cHCName[iDet]);
115 // histograms IDs
116 fDriftHistoID[kH1][iDet] = analysisManager->GetH1Id(histoName[kH1][iDet]);
117 fDriftHistoID[kH2][iDet] = analysisManager->GetH2Id(histoName[kH2][iDet]);
118 }
119 }
120}
121
122//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
123
125{
126 //
127 // Fill histograms & ntuple
128 //
129
130 // Get analysis manager
131 auto analysisManager = G4AnalysisManager::Instance();
132
133 // Drift chambers hits
134 for (G4int iDet = 0; iDet < kDim; ++iDet) {
135 auto hc = GetHC(event, fDriftHCID[iDet]);
136 if ( ! hc ) return;
137
138 auto nhit = hc->GetSize();
139 analysisManager->FillH1(fDriftHistoID[kH1][iDet], nhit );
140 // columns 0, 1
141 analysisManager->FillNtupleIColumn(iDet, nhit);
142
143 for (unsigned long i = 0; i < nhit; ++i) {
144 auto hit = static_cast<DriftChamberHit*>(hc->GetHit(i));
145 auto localPos = hit->GetLocalPos();
146 analysisManager->FillH2(fDriftHistoID[kH2][iDet], localPos.x(), localPos.y());
147 }
148 }
149
150 // Em/Had Calorimeters hits
151 array<G4int, kDim> totalCalHit = {{ 0, 0 }};
152 array<G4double, kDim> totalCalEdep = {{ 0., 0. }};
153
154 for (G4int iDet = 0; iDet < kDim; ++iDet) {
155 auto hc = GetHC(event, fCalHCID[iDet]);
156 if ( ! hc ) return;
157
158 totalCalHit[iDet] = 0;
159 totalCalEdep[iDet] = 0.;
160 for (unsigned long i = 0; i < hc->GetSize(); ++i) {
161 G4double edep = 0.;
162 // The EM and Had calorimeter hits are of different types
163 if (iDet == 0) {
164 auto hit = static_cast<EmCalorimeterHit*>(hc->GetHit(i));
165 edep = hit->GetEdep();
166 } else {
167 auto hit = static_cast<HadCalorimeterHit*>(hc->GetHit(i));
168 edep = hit->GetEdep();
169 }
170 if ( edep > 0. ) {
171 totalCalHit[iDet]++;
172 totalCalEdep[iDet] += edep;
173 }
174 fCalEdep[iDet][i] = edep;
175 }
176 // columns 2, 3
177 analysisManager->FillNtupleDColumn(iDet + 2, totalCalEdep[iDet]);
178 }
179
180 // Hodoscopes hits
181 for (G4int iDet = 0; iDet < kDim; ++iDet) {
182 auto hc = GetHC(event, fHodHCID[iDet]);
183 if ( ! hc ) return;
184
185 for (unsigned int i = 0; i<hc->GetSize(); ++i) {
186 auto hit = static_cast<HodoscopeHit*>(hc->GetHit(i));
187 // columns 4, 5
188 analysisManager->FillNtupleDColumn(iDet + 4, hit->GetTime());
189 }
190 }
191 analysisManager->AddNtupleRow();
192
193 //
194 // Print diagnostics
195 //
196
197 auto printModulo = G4RunManager::GetRunManager()->GetPrintProgress();
198 if ( printModulo == 0 || event->GetEventID() % printModulo != 0) return;
199
200 auto primary = event->GetPrimaryVertex(0)->GetPrimary(0);
201 G4cout
202 << G4endl
203 << ">>> Event " << event->GetEventID() << " >>> Simulation truth : "
204 << primary->GetG4code()->GetParticleName()
205 << " " << primary->GetMomentum() << G4endl;
206
207 // Hodoscopes
208 for (G4int iDet = 0; iDet < kDim; ++iDet) {
209 auto hc = GetHC(event, fHodHCID[iDet]);
210 if ( ! hc ) return;
211 G4cout << "Hodoscope " << iDet + 1 << " has " << hc->GetSize() << " hits." << G4endl;
212 for (unsigned int i = 0; i<hc->GetSize(); ++i) {
213 hc->GetHit(i)->Print();
214 }
215 }
216
217 // Drift chambers
218 for (G4int iDet = 0; iDet < kDim; ++iDet) {
219 auto hc = GetHC(event, fDriftHCID[iDet]);
220 if ( ! hc ) return;
221 G4cout << "Drift Chamber " << iDet + 1 << " has " << hc->GetSize() << " hits." << G4endl;
222 for (auto layer = 0; layer < kNofChambers; ++layer) {
223 for (unsigned int i = 0; i < hc->GetSize(); i++) {
224 auto hit = static_cast<DriftChamberHit*>(hc->GetHit(i));
225 if (hit->GetLayerID() == layer) hit->Print();
226 }
227 }
228 }
229
230 // Calorimeters
231 array<G4String, kDim> calName = {{ "EM", "Hadron" }};
232 for (G4int iDet = 0; iDet < kDim; ++iDet) {
233 G4cout << calName[iDet] << " Calorimeter has " << totalCalHit[iDet] << " hits."
234 << " Total Edep is " << totalCalEdep[iDet]/MeV << " (MeV)" << G4endl;
235 }
236}
237
238//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
239
240}
Definition of B5 example constants.
Definition of the B5::DriftChamberHit class.
Definition of the B5::EmCalorimeterHit class.
Definition of the B5::HadCalorimeterHit class.
Definition of the B5::HodoscopeHit class.
const G4int kH1
const G4int kH2
const G4int kDim
Drift chamber hit.
G4ThreeVector GetLocalPos() const
void Print() override
EM Calorimeter hit.
G4double GetEdep() const
std::array< G4int, kDim > fDriftHCID
std::array< std::vector< G4double >, kDim > fCalEdep
void EndOfEventAction(const G4Event *) override
std::array< std::array< G4int, kDim >, kDim > fDriftHistoID
std::array< G4int, kDim > fCalHCID
std::array< G4int, kDim > fHodHCID
void BeginOfEventAction(const G4Event *) override
Hadron Calorimeter hit.
constexpr G4int kNofChambers
Definition Constants.hh:40

Applications | User Support | Publications | Collaboration