Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
B5EventAction.cc
이 파일의 문서화 페이지로 가기
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: B5EventAction.cc 108609 2018-02-26 10:22:22Z gcosmo $
27 //
30 
31 #include "B5EventAction.hh"
32 #include "B5HodoscopeHit.hh"
33 #include "B5DriftChamberHit.hh"
34 #include "B5EmCalorimeterHit.hh"
35 #include "B5HadCalorimeterHit.hh"
36 #include "B5Constants.hh"
37 #include "B5Analysis.hh"
38 
39 #include "G4Event.hh"
40 #include "G4RunManager.hh"
41 #include "G4EventManager.hh"
42 #include "G4HCofThisEvent.hh"
43 #include "G4VHitsCollection.hh"
44 #include "G4SDManager.hh"
45 #include "G4SystemOfUnits.hh"
46 #include "G4ios.hh"
47 
48 using std::array;
49 using std::vector;
50 
51 
52 namespace {
53 
54 // Utility function which finds a hit collection with the given Id
55 // and print warnings if not found
56 G4VHitsCollection* GetHC(const G4Event* event, G4int collId) {
57  auto hce = event->GetHCofThisEvent();
58  if (!hce) {
60  msg << "No hits collection of this event found." << G4endl;
61  G4Exception("B5EventAction::EndOfEventAction()",
62  "B5Code001", JustWarning, msg);
63  return nullptr;
64  }
65 
66  auto hc = hce->GetHC(collId);
67  if ( ! hc) {
69  msg << "Hits collection " << collId << " of this event not found." << G4endl;
70  G4Exception("B5EventAction::EndOfEventAction()",
71  "B5Code001", JustWarning, msg);
72  }
73  return hc;
74 }
75 
76 }
77 
78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
79 
82  fHodHCID {{ -1, -1 }},
83  fDriftHCID{{ -1, -1 }},
84  fCalHCID {{ -1, -1 }},
85  fDriftHistoID{{ {{ -1, -1 }}, {{ -1, -1 }} }},
86  fCalEdep{{ vector<G4double>(kNofEmCells, 0.), vector<G4double>(kNofHadCells, 0.) }}
87  // std::array<T, N> is an aggregate that contains a C array.
88  // To initialize it, we need outer braces for the class itself
89  // and inner braces for the C array
90 {
91  // set printing per each event
93 }
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 
98 {}
99 
100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
101 
103 {
104  // Find hit collections and histogram Ids by names (just once)
105  // and save them in the data members of this class
106 
107  if (fHodHCID[0] == -1) {
108  auto sdManager = G4SDManager::GetSDMpointer();
109  auto analysisManager = G4AnalysisManager::Instance();
110 
111  // hits collections names
112  array<G4String, kDim> hHCName
113  = {{ "hodoscope1/hodoscopeColl", "hodoscope2/hodoscopeColl" }};
114  array<G4String, kDim> dHCName
115  = {{ "chamber1/driftChamberColl", "chamber2/driftChamberColl" }};
116  array<G4String, kDim> cHCName
117  = {{ "EMcalorimeter/EMcalorimeterColl", "HadCalorimeter/HadCalorimeterColl" }};
118 
119  // histograms names
120  array<array<G4String, kDim>, kDim> histoName
121  = {{ {{ "Chamber1", "Chamber2" }}, {{ "Chamber1 XY", "Chamber2 XY" }} }};
122 
123  for (G4int iDet = 0; iDet < kDim; ++iDet) {
124  // hit collections IDs
125  fHodHCID[iDet] = sdManager->GetCollectionID(hHCName[iDet]);
126  fDriftHCID[iDet] = sdManager->GetCollectionID(dHCName[iDet]);
127  fCalHCID[iDet] = sdManager->GetCollectionID(cHCName[iDet]);
128  // histograms IDs
129  fDriftHistoID[kH1][iDet] = analysisManager->GetH1Id(histoName[kH1][iDet]);
130  fDriftHistoID[kH2][iDet] = analysisManager->GetH2Id(histoName[kH2][iDet]);
131  }
132  }
133 }
134 
135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136 
137 void B5EventAction::EndOfEventAction(const G4Event* event)
138 {
139  //
140  // Fill histograms & ntuple
141  //
142 
143  // Get analysis manager
144  auto analysisManager = G4AnalysisManager::Instance();
145 
146  // Drift chambers hits
147  for (G4int iDet = 0; iDet < kDim; ++iDet) {
148  auto hc = GetHC(event, fDriftHCID[iDet]);
149  if ( ! hc ) return;
150 
151  auto nhit = hc->GetSize();
152  analysisManager->FillH1(fDriftHistoID[kH1][iDet], nhit );
153  // columns 0, 1
154  analysisManager->FillNtupleIColumn(iDet, nhit);
155 
156  for (unsigned long i = 0; i < nhit; ++i) {
157  auto hit = static_cast<B5DriftChamberHit*>(hc->GetHit(i));
158  auto localPos = hit->GetLocalPos();
159  analysisManager->FillH2(fDriftHistoID[kH2][iDet], localPos.x(), localPos.y());
160  }
161  }
162 
163  // Em/Had Calorimeters hits
164  array<G4int, kDim> totalCalHit = {{ 0, 0 }};
165  array<G4double, kDim> totalCalEdep = {{ 0., 0. }};
166 
167  for (G4int iDet = 0; iDet < kDim; ++iDet) {
168  auto hc = GetHC(event, fCalHCID[iDet]);
169  if ( ! hc ) return;
170 
171  totalCalHit[iDet] = 0;
172  totalCalEdep[iDet] = 0.;
173  for (unsigned long i = 0; i < hc->GetSize(); ++i) {
174  G4double edep = 0.;
175  // The EM and Had calorimeter hits are of different types
176  if (iDet == 0) {
177  auto hit = static_cast<B5EmCalorimeterHit*>(hc->GetHit(i));
178  edep = hit->GetEdep();
179  } else {
180  auto hit = static_cast<B5HadCalorimeterHit*>(hc->GetHit(i));
181  edep = hit->GetEdep();
182  }
183  if ( edep > 0. ) {
184  totalCalHit[iDet]++;
185  totalCalEdep[iDet] += edep;
186  }
187  fCalEdep[iDet][i] = edep;
188  }
189  // columns 2, 3
190  analysisManager->FillNtupleDColumn(iDet + 2, totalCalEdep[iDet]);
191  }
192 
193  // Hodoscopes hits
194  for (G4int iDet = 0; iDet < kDim; ++iDet) {
195  auto hc = GetHC(event, fHodHCID[iDet]);
196  if ( ! hc ) return;
197 
198  for (unsigned int i = 0; i<hc->GetSize(); ++i) {
199  auto hit = static_cast<B5HodoscopeHit*>(hc->GetHit(i));
200  // columns 4, 5
201  analysisManager->FillNtupleDColumn(iDet + 4, hit->GetTime());
202  }
203  }
204  analysisManager->AddNtupleRow();
205 
206  //
207  // Print diagnostics
208  //
209 
210  auto printModulo = G4RunManager::GetRunManager()->GetPrintProgress();
211  if ( printModulo == 0 || event->GetEventID() % printModulo != 0) return;
212 
213  auto primary = event->GetPrimaryVertex(0)->GetPrimary(0);
214  G4cout
215  << G4endl
216  << ">>> Event " << event->GetEventID() << " >>> Simulation truth : "
217  << primary->GetG4code()->GetParticleName()
218  << " " << primary->GetMomentum() << G4endl;
219 
220  // Hodoscopes
221  for (G4int iDet = 0; iDet < kDim; ++iDet) {
222  auto hc = GetHC(event, fHodHCID[iDet]);
223  if ( ! hc ) return;
224  G4cout << "Hodoscope " << iDet + 1 << " has " << hc->GetSize() << " hits." << G4endl;
225  for (unsigned int i = 0; i<hc->GetSize(); ++i) {
226  hc->GetHit(i)->Print();
227  }
228  }
229 
230  // Drift chambers
231  for (G4int iDet = 0; iDet < kDim; ++iDet) {
232  auto hc = GetHC(event, fDriftHCID[iDet]);
233  if ( ! hc ) return;
234  G4cout << "Drift Chamber " << iDet + 1 << " has " << hc->GetSize() << " hits." << G4endl;
235  for (auto layer = 0; layer < kNofChambers; ++layer) {
236  for (unsigned int i = 0; i < hc->GetSize(); i++) {
237  auto hit = static_cast<B5DriftChamberHit*>(hc->GetHit(i));
238  if (hit->GetLayerID() == layer) hit->Print();
239  }
240  }
241  }
242 
243  // Calorimeters
244  array<G4String, kDim> calName = {{ "EM", "Hadron" }};
245  for (G4int iDet = 0; iDet < kDim; ++iDet) {
246  G4cout << calName[iDet] << " Calorimeter has " << totalCalHit[iDet] << " hits."
247  << " Total Edep is " << totalCalEdep[iDet]/MeV << " (MeV)" << G4endl;
248  }
249 }
250 
251 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Definition of the B5DriftChamberHit class.
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:80
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:45
static constexpr double MeV
Definition: G4SIunits.hh:214
G4int GetEventID() const
Definition: G4Event.hh:151
#define G4endl
Definition: G4ios.hh:61
Definition of the B5HadCalorimeterHit class.
const G4int kH1
std::array< G4int, kDim > fHodHCID
Selection of the analysis technology.
fCalHCID
G4double GetEdep() const
virtual ~B5EventAction()
constexpr G4int kNofHadCells
Definition: B5Constants.hh:44
double G4double
Definition: G4Types.hh:76
Definition of the B5EventAction class.
const G4int kDim
Definition: B4bRunData.hh:43
fDriftHCID
Double_t edep
std::array< std::vector< G4double >, kDim > fCalEdep
G4int GetPrintProgress()
G4double GetEdep() const
std::array< std::array< G4int, kDim >, kDim > fDriftHistoID
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.hh:65
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
std::array< G4int, kDim > fDriftHCID
int G4int
Definition: G4Types.hh:78
std::array< G4int, kDim > fCalHCID
const G4int kH2
constexpr G4int kNofChambers
Definition: B5Constants.hh:38
fDriftHistoID
G4GLOB_DLL std::ostream G4cout
virtual void Print()
Definition of the B5HodoscopeHit class.
virtual void BeginOfEventAction(const G4Event *)
virtual void EndOfEventAction(const G4Event *)
Definition of constants.
constexpr G4int kNofEmCells
Definition: B5Constants.hh:41
const G4double hc
[MeV*fm]
Definition of the B5EmCalorimeterHit class.
G4ThreeVector GetLocalPos() const
void SetPrintProgress(G4int i)