Loading...
Searching...
No Matches
Par04EventAction.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#include "Par04EventAction.hh"
27
28#include "Par04DetectorConstruction.hh" // for Par04DetectorConstruction
29#include "Par04Hit.hh" // for Par04Hit, Par04HitsCollection
31
32#include "G4AnalysisManager.hh" // for G4AnalysisManager
33#include "G4Event.hh" // for G4Event
34#include "G4EventManager.hh" // for G4EventManager
35#include "G4HCofThisEvent.hh" // for G4HCofThisEvent
36#include "G4SDManager.hh" // for G4SDManager
37
38#include <CLHEP/Units/SystemOfUnits.h> // for GeV
39#include <CLHEP/Vector/ThreeVector.h> // for Hep3Vector
40#include "G4Exception.hh" // for G4Exception, G4ExceptionDesc...
41#include "G4ExceptionSeverity.hh" // for FatalException
42#include "G4GenericAnalysisManager.hh" // for G4GenericAnalysisManager
43#include "G4PrimaryParticle.hh" // for G4PrimaryParticle
44#include "G4PrimaryVertex.hh" // for G4PrimaryVertex
45#include "G4SystemOfUnits.hh" // for GeV
46#include "G4THitsCollection.hh" // for G4THitsCollection
47#include "G4ThreeVector.hh" // for G4ThreeVector
48#include "G4Timer.hh" // for G4Timer
49#include "G4UserEventAction.hh" // for G4UserEventAction
50#include <algorithm> // for max
51#include <cmath> // for log10
52#include <cstddef> // for size_t
53#include <ostream> // for basic_ostream::operator<<
54
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
56
58 Par04ParallelFullWorld* aParallel)
60 , fHitCollectionID(-1)
61 , fPhysicalFullHitCollectionID(-1)
62 , fPhysicalFastHitCollectionID(-1)
63 , fTimer()
64 , fDetector(aDetector)
65 , fParallel(aParallel)
66{
67 fCellNbRho = aDetector->GetMeshNbOfCells().x();
68 fCellNbPhi = aDetector->GetMeshNbOfCells().y();
69 fCellNbZ = aDetector->GetMeshNbOfCells().z();
78}
79
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
81
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85
89
90//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
91
93 fTimer.Start();
94}
95
96//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
97
99 fTimer.Stop();
100}
101
102//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103
105{
106 G4SDManager::GetSDMpointer()->GetHCtable();
107 StopTimer();
108
109 // Get hits collection ID (only once)
110 if(fHitCollectionID == -1)
111 {
112 fHitCollectionID = G4SDManager::GetSDMpointer()->GetCollectionID("hits");
113 }
115 {
117 G4SDManager::GetSDMpointer()->GetCollectionID("physicalCellsFullSim");
118 }
120 {
122 G4SDManager::GetSDMpointer()->GetCollectionID("physicalCellsFastSim");
123 }
124 // Get hits collection
125 auto hitsCollection =
126 static_cast<Par04HitsCollection*>(aEvent->GetHCofThisEvent()->GetHC(fHitCollectionID));
127 auto physicalFullHitsCollection =
128 static_cast<Par04HitsCollection*>(aEvent->GetHCofThisEvent()
130 auto physicalFastHitsCollection =
131 static_cast<Par04HitsCollection*>(aEvent->GetHCofThisEvent()
133
134 if(hitsCollection == nullptr)
135 {
136 G4ExceptionDescription msg;
137 msg << "Cannot access hitsCollection ID " << fHitCollectionID;
138 G4Exception("Par04EventAction::GetHitsCollection()", "MyCode0001", FatalException, msg);
139 }
140 if(physicalFullHitsCollection == nullptr)
141 {
142 G4ExceptionDescription msg;
143 msg << "Cannot access physical full sim hitsCollection ID " << fPhysicalFullHitCollectionID;
144 G4Exception("Par04EventAction::GetHitsCollection()", "MyCode0001", FatalException, msg);
145 }
146 if(physicalFastHitsCollection == nullptr)
147 {
148 G4ExceptionDescription msg;
149 msg << "Cannot access physical fast sim hitsCollection ID " << fPhysicalFastHitCollectionID;
150 G4Exception("Par04EventAction::GetHitsCollection()", "MyCode0001", FatalException, msg);
151 }
152 // Get analysis manager
153 auto analysisManager = G4AnalysisManager::Instance();
154 // Retrieve only once detector dimensions
155 if(fCellSizeZ == 0)
156 {
163 }
164 if(fPhysicalNbLayers == 0)
165 {
169 }
170
171 // Retrieve information from primary vertex and primary particle
172 // To calculate shower axis and entry point to the detector
173 auto primaryVertex =
174 G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetPrimaryVertex();
175 auto primaryParticle = primaryVertex->GetPrimary(0);
176 G4double primaryEnergy = primaryParticle->GetTotalEnergy();
177 // Estimate from vertex and particle direction the entry point to the detector
178 // Calculate entrance point to the detector located at z = 0
179 auto primaryDirection = primaryParticle->GetMomentumDirection();
180 auto primaryEntrance =
181 primaryVertex->GetPosition() - primaryVertex->GetPosition().z() * primaryDirection;
182
183 // Resize back to initial mesh size
187 fCalZ.resize(fCellNbRho * fCellNbPhi * fCellNbZ,0);
192
193 // Fill histograms
194 Par04Hit* hit = nullptr;
195 G4double hitEn = 0;
196 G4double totalEnergy = 0;
197 G4int hitNum = 0;
198 G4int totalNum = 0;
199 G4int hitZ = -1;
200 G4int hitRho = -1;
201 G4int hitPhi = -1;
202 G4int hitType = -1;
203 G4int numNonZeroThresholdCells = 0;
204 G4double tDistance = 0., rDistance = 0., phiDistance = 0.;
205 G4double tFirstMoment = 0., tSecondMoment = 0.;
206 G4double rFirstMoment = 0., rSecondMoment = 0.;
207 G4double phiMean = 0.;
208 for(size_t iHit = 0; iHit < hitsCollection->entries(); iHit++)
209 {
210 hit = static_cast<Par04Hit*>(hitsCollection->GetHit(iHit));
211 hitZ = hit->GetZid();
212 hitRho = hit->GetRhoId();
213 hitPhi = hit->GetPhiId();
214 hitEn = hit->GetEdep();
215 hitNum = hit->GetNdep();
216 hitType = hit->GetType();
217 if(hitEn > 0)
218 {
219 totalEnergy += hitEn;
220 totalNum += hitNum;
221 tDistance = hitZ * fCellSizeZ;
222 rDistance = hitRho * fCellSizeRho;
223 phiDistance = hitPhi * fCellSizePhi;
224 tFirstMoment += hitEn * tDistance;
225 rFirstMoment += hitEn * rDistance;
226 phiMean += hitEn * phiDistance;
227 analysisManager->FillH1(4, tDistance, hitEn);
228 analysisManager->FillH1(5, rDistance, hitEn);
229 analysisManager->FillH1(10, hitType);
230 if(hitEn > 0.0005)
231 { // e > 0.5 keV
232 fCalEdep[numNonZeroThresholdCells] = hitEn;
233 fCalRho[numNonZeroThresholdCells] = hitRho;
234 fCalPhi[numNonZeroThresholdCells] = hitPhi;
235 fCalZ[numNonZeroThresholdCells] = hitZ;
236 numNonZeroThresholdCells++;
237 analysisManager->FillH1(13, std::log10(hitEn));
238 analysisManager->FillH1(15, hitNum);
239 }
240 }
241 }
242 tFirstMoment /= totalEnergy;
243 rFirstMoment /= totalEnergy;
244 phiMean /= totalEnergy;
245 analysisManager->FillH1(0, primaryEnergy / GeV);
246 analysisManager->FillH1(1, totalEnergy / GeV);
247 analysisManager->FillH1(2, totalEnergy / primaryEnergy);
248 analysisManager->FillH1(3, fTimer.GetRealElapsed());
249 analysisManager->FillH1(6, tFirstMoment);
250 analysisManager->FillH1(7, rFirstMoment);
251 analysisManager->FillH1(12, numNonZeroThresholdCells);
252 analysisManager->FillH1(14, totalNum);
253 // Resize to store only energy hits above threshold
254 fCalEdep.resize(numNonZeroThresholdCells);
255 fCalRho.resize(numNonZeroThresholdCells);
256 fCalPhi.resize(numNonZeroThresholdCells);
257 fCalZ.resize(numNonZeroThresholdCells);
258 analysisManager->FillNtupleDColumn(0, 0, primaryEnergy);
259 analysisManager->FillNtupleDColumn(0, 1, fTimer.GetRealElapsed());
260 // Second loop over hits to calculate second moments
261 for(size_t iHit = 0; iHit < hitsCollection->entries(); iHit++)
262 {
263 hit = static_cast<Par04Hit*>(hitsCollection->GetHit(iHit));
264 hitEn = hit->GetEdep();
265 hitZ = hit->GetZid();
266 hitRho = hit->GetRhoId();
267 hitPhi = hit->GetPhiId();
268 if(hitEn > 0)
269 {
270 tDistance = hitZ * fCellSizeZ;
271 rDistance = hitRho * fCellSizeRho;
272 phiDistance = hitPhi * fCellSizePhi;
273 tSecondMoment += hitEn * std::pow(tDistance - tFirstMoment, 2);
274 rSecondMoment += hitEn * std::pow(rDistance - rFirstMoment, 2);
275 analysisManager->FillH1(11, phiDistance - phiMean, hitEn);
276 }
277 }
278 tSecondMoment /= totalEnergy;
279 rSecondMoment /= totalEnergy;
280 analysisManager->FillH1(8, tSecondMoment);
281 analysisManager->FillH1(9, rSecondMoment);
282
283 // Fill ntuple with physical readout data
284 G4double totalPhysicalEnergy = 0;
285 totalNum = 0;
286 hitEn = 0;
287 hitNum = 0;
288 G4int hitLayer = -1;
289 G4int hitRow = -1;
290 G4int hitSlice = -1;
291 numNonZeroThresholdCells = 0;
292 for(size_t iHit = 0; iHit < physicalFullHitsCollection->entries(); iHit++)
293 {
294 hit = static_cast<Par04Hit*>(physicalFullHitsCollection->GetHit(iHit));
295 hitLayer = hit->GetRhoId();
296 hitRow = hit->GetZid();
297 hitSlice = hit->GetPhiId();
298 hitEn = hit->GetEdep();
299 hitNum = hit->GetNdep();
300 if(hitEn > 0)
301 {
302 totalPhysicalEnergy += hitEn;
303 totalNum += hitNum;
304 if(hitEn > 0.0005)
305 { // e > 0.5 keV
306 fCalPhysicalEdep[numNonZeroThresholdCells] = hitEn;
307 fCalPhysicalLayer[numNonZeroThresholdCells] = hitLayer;
308 fCalPhysicalRow[numNonZeroThresholdCells] = hitRow;
309 fCalPhysicalSlice[numNonZeroThresholdCells] = hitSlice;
310 numNonZeroThresholdCells++;
311 analysisManager->FillH1(19, std::log10(hitEn));
312 analysisManager->FillH1(21, hitNum);
313 }
314 }
315 }for(size_t iHit = 0; iHit < physicalFastHitsCollection->entries(); iHit++)
316 {
317 hit = static_cast<Par04Hit*>(physicalFastHitsCollection->GetHit(iHit));
318 hitLayer = hit->GetRhoId();
319 hitRow = hit->GetZid();
320 hitSlice = hit->GetPhiId();
321 hitEn = hit->GetEdep();
322 hitNum = hit->GetNdep();
323 if(hitEn > 0)
324 {
325 totalPhysicalEnergy += hitEn;
326 totalNum += hitNum;
327 if(hitEn > 0.0005)
328 { // e > 0.5 keV
329 fCalPhysicalEdep[numNonZeroThresholdCells] = hitEn;
330 fCalPhysicalLayer[numNonZeroThresholdCells] = hitLayer;
331 fCalPhysicalRow[numNonZeroThresholdCells] = hitRow;
332 fCalPhysicalSlice[numNonZeroThresholdCells] = hitSlice;
333 numNonZeroThresholdCells++;
334 analysisManager->FillH1(19, std::log10(hitEn));
335 analysisManager->FillH1(21, hitNum);
336 }
337 }
338 }
339 analysisManager->FillH1(16, totalPhysicalEnergy / GeV);
340 analysisManager->FillH1(17, totalPhysicalEnergy / primaryEnergy);
341 analysisManager->FillH1(18, numNonZeroThresholdCells);
342 analysisManager->FillH1(20, totalNum);
343 fCalPhysicalEdep.resize(numNonZeroThresholdCells);
344 fCalPhysicalLayer.resize(numNonZeroThresholdCells);
345 fCalPhysicalSlice.resize(numNonZeroThresholdCells);
346 fCalPhysicalRow.resize(numNonZeroThresholdCells);
347 analysisManager->AddNtupleRow(0);
348 analysisManager->AddNtupleRow(1);
349 analysisManager->AddNtupleRow(2);
350}
G4THitsCollection< Par04Hit > Par04HitsCollection
Definition Par04Hit.hh:145
G4ThreeVector GetMeshSizeOfCells() const
Get size of Mesh cells in cylindrical coordinates (r, phi, z)
G4ThreeVector GetMeshNbOfCells() const
Get number of Mesh cells in cylindrical coordinates (r, phi, z)
std::vector< G4int > fCalZ
Cell ID of z axis to be stored in ntuple.
G4double fCellSizePhi
Size of cell in azimuthal angle.
virtual void BeginOfEventAction(const G4Event *aEvent) final
Timer is started.
G4int fHitCollectionID
ID of a hit collection to analyse.
std::vector< G4int > fCalPhi
Cell ID of azimuthal angle to be stored in ntuple.
G4int fCellNbRho
Number of readout cells along radius.
std::vector< G4int > fCalPhysicalRow
Physical row ID to be stored in ntuple.
virtual ~Par04EventAction()
G4int fCellNbZ
Number of readout cells along z axis.
G4int fPhysicalNbRows
Number of physical readout rows.
std::vector< G4double > fCalPhysicalEdep
Physical cell energy deposits to be stored in ntuple.
G4double fCellSizeRho
Size of cell along radius of cylinder.
Par04DetectorConstruction * fDetector
Pointer to detector construction to retrieve (once) the detector dimensions and size of readout.
G4int fPhysicalNbLayers
Number of physical readout layers.
G4double fCellSizeZ
Size of cell along Z axis.
std::vector< G4int > fCalRho
Cell ID of radius to be stored in ntuple.
G4Timer fTimer
Timer measurement from Geant4.
G4int fPhysicalNbSlices
Number of physical readout slices.
std::vector< G4double > fCalEdep
Cell energy deposits to be stored in ntuple.
Par04EventAction(Par04DetectorConstruction *aDetector, Par04ParallelFullWorld *aParallel)
std::vector< G4int > fCalPhysicalLayer
Physical layer ID to be stored in ntuple.
virtual void EndOfEventAction(const G4Event *aEvent) final
Hits collection is retrieved, analysed, and histograms are filled.
std::vector< G4int > fCalPhysicalSlice
Physical slice ID to be stored in ntuple.
G4int fCellNbPhi
Number of readout cells in azimuthal angle.
Par04ParallelFullWorld * fParallel
Hit class to store energy deposited in the sensitive detector.
G4int GetZid() const
Get Z id of the cell in the readout segmentation.
Definition Par04Hit.hh:100
G4double GetEdep() const
Get energy.
Definition Par04Hit.hh:90
G4int GetType() const
Get type (0 = full sim, 1 = fast sim)
Definition Par04Hit.hh:116
G4int GetRhoId() const
Get rho id of the cell in the readout segmentation.
Definition Par04Hit.hh:104
G4int GetNdep() const
Get number of deposits per hit/cell.
Definition Par04Hit.hh:96
G4int GetPhiId() const
Get phi id of the cell in the readout segmentation.
Definition Par04Hit.hh:108
G4int GetNbOfRows() const
Get number of rows.
G4int GetNbOfLayers() const
Get number of layers.
G4int GetNbOfSlices() const
Get number of slices.

Applications | User Support | Publications | Collaboration