Loading...
Searching...
No Matches
LXeEventAction.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 optical/LXe/src/LXeEventAction.cc
28/// \brief Implementation of the LXeEventAction class
29//
30//
31#include "LXeEventAction.hh"
32
34#include "LXeHistoManager.hh"
35#include "LXePMTHit.hh"
36#include "LXeRun.hh"
37#include "LXeScintHit.hh"
38#include "LXeTrajectory.hh"
39
40#include "G4Event.hh"
41#include "G4EventManager.hh"
42#include "G4ios.hh"
43#include "G4RunManager.hh"
44#include "G4SDManager.hh"
45#include "G4SystemOfUnits.hh"
46#include "G4Trajectory.hh"
47#include "G4TrajectoryContainer.hh"
48#include "G4UImanager.hh"
49#include "G4VVisManager.hh"
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52
54 : fDetector(det)
55{
57}
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
60
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
64
66{
67 fHitCount = 0;
72 fTotE = 0.0;
73
74 fConvPosSet = false;
75 fEdepMax = 0.0;
76
78
79 G4SDManager* SDman = G4SDManager::GetSDMpointer();
80 if(fScintCollID < 0)
81 fScintCollID = SDman->GetCollectionID("scintCollection");
82 if(fPMTCollID < 0)
83 fPMTCollID = SDman->GetCollectionID("pmtHitCollection");
84}
85
86//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87
89{
90 G4TrajectoryContainer* trajectoryContainer =
91 anEvent->GetTrajectoryContainer();
92
93 G4int n_trajectories = 0;
94 if(trajectoryContainer)
95 n_trajectories = trajectoryContainer->entries();
96
97 // extract the trajectories and draw them
98 if(G4VVisManager::GetConcreteInstance())
99 {
100 for(G4int i = 0; i < n_trajectories; ++i)
101 {
102 auto trj = (LXeTrajectory*) ((*(anEvent->GetTrajectoryContainer()))[i]);
103 if(trj->GetParticleName() == "opticalphoton")
104 {
105 trj->SetForceDrawTrajectory(fForcedrawphotons);
106 trj->SetForceNoDrawTrajectory(fForcenophotons);
107 }
108 trj->DrawTrajectory();
109 }
110 }
111
112 LXeScintHitsCollection* scintHC = nullptr;
113 LXePMTHitsCollection* pmtHC = nullptr;
114 G4HCofThisEvent* hitsCE = anEvent->GetHCofThisEvent();
115
116 // Get the hit collections
117 if(hitsCE)
118 {
119 if(fScintCollID >= 0)
120 {
121 scintHC = (LXeScintHitsCollection*) (hitsCE->GetHC(fScintCollID));
122 }
123 if(fPMTCollID >= 0)
124 {
125 pmtHC = (LXePMTHitsCollection*) (hitsCE->GetHC(fPMTCollID));
126 }
127 }
128
129 // Hits in scintillator
130 if(scintHC)
131 {
132 size_t n_hit = scintHC->entries();
133 G4ThreeVector eWeightPos(0.);
134 G4double edep;
135 G4double edepMax = 0;
136
137 for(size_t i = 0; i < n_hit; ++i)
138 { // gather info on hits in scintillator
139 edep = (*scintHC)[i]->GetEdep();
140 fTotE += edep;
141 eWeightPos +=
142 (*scintHC)[i]->GetPos() * edep; // calculate energy weighted pos
143 if(edep > edepMax)
144 {
145 edepMax = edep; // store max energy deposit
146 G4ThreeVector posMax = (*scintHC)[i]->GetPos();
147 fPosMax = posMax;
148 fEdepMax = edep;
149 }
150 }
151
152 G4AnalysisManager::Instance()->FillH1(7, fTotE);
153
154 if(fTotE == 0.)
155 {
156 if(fVerbose > 0)
157 G4cout << "No hits in the scintillator this event." << G4endl;
158 }
159 else
160 {
161 // Finish calculation of energy weighted position
162 eWeightPos /= fTotE;
163 fEWeightPos = eWeightPos;
164 if(fVerbose > 0)
165 {
166 G4cout << "\tEnergy weighted position of hits in LXe : "
167 << eWeightPos / mm << G4endl;
168 }
169 }
170 if(fVerbose > 0)
171 {
172 G4cout << "\tTotal energy deposition in scintillator : " << fTotE / keV
173 << " (keV)" << G4endl;
174 }
175 }
176
177 if(pmtHC)
178 {
179 G4ThreeVector reconPos(0., 0., 0.);
180 size_t pmts = pmtHC->entries();
181 // Gather info from all PMTs
182 for(size_t i = 0; i < pmts; ++i)
183 {
184 fHitCount += (*pmtHC)[i]->GetPhotonCount();
185 reconPos += (*pmtHC)[i]->GetPMTPos() * (*pmtHC)[i]->GetPhotonCount();
186 if((*pmtHC)[i]->GetPhotonCount() >= fPMTThreshold)
187 {
189 }
190 else
191 { // wasn't above the threshold, turn it back off
192 (*pmtHC)[i]->SetDrawit(false);
193 }
194 }
195
196 G4AnalysisManager::Instance()->FillH1(1, fHitCount);
197 G4AnalysisManager::Instance()->FillH1(2, fPMTsAboveThreshold);
198
199 if(fHitCount > 0)
200 { // don't bother unless there were hits
201 reconPos /= fHitCount;
202 if(fVerbose > 0)
203 {
204 G4cout << "\tReconstructed position of hits in LXe : " << reconPos / mm
205 << G4endl;
206 }
207 fReconPos = reconPos;
208 }
209 pmtHC->DrawAllHits();
210 }
211
212 G4AnalysisManager::Instance()->FillH1(3, fPhotonCount_Scint);
213 G4AnalysisManager::Instance()->FillH1(4, fPhotonCount_Ceren);
214 G4AnalysisManager::Instance()->FillH1(5, fAbsorptionCount);
215 G4AnalysisManager::Instance()->FillH1(6, fBoundaryAbsorptionCount);
216
217 if(fVerbose > 0)
218 {
219 // End of event output. later to be controlled by a verbose level
220 G4cout << "\tNumber of photons that hit PMTs in this event : " << fHitCount
221 << G4endl;
222 G4cout << "\tNumber of PMTs above threshold(" << fPMTThreshold
223 << ") : " << fPMTsAboveThreshold << G4endl;
224 G4cout << "\tNumber of photons produced by scintillation in this event : "
225 << fPhotonCount_Scint << G4endl;
226 G4cout << "\tNumber of photons produced by cerenkov in this event : "
227 << fPhotonCount_Ceren << G4endl;
228 G4cout << "\tNumber of photons absorbed (OpAbsorption) in this event : "
229 << fAbsorptionCount << G4endl;
230 G4cout << "\tNumber of photons absorbed at boundaries (OpBoundary) in "
231 << "this event : " << fBoundaryAbsorptionCount << G4endl;
232 G4cout << "Unaccounted for photons in this event : "
235 << G4endl;
236 }
237
238 // update the run statistics
239 auto run = static_cast<LXeRun*>(
240 G4RunManager::GetRunManager()->GetNonConstCurrentRun());
241
243 run->IncPhotonCount_Scint(fPhotonCount_Scint);
244 run->IncPhotonCount_Ceren(fPhotonCount_Ceren);
245 run->IncEDep(fTotE);
246 run->IncAbsorption(fAbsorptionCount);
247 run->IncBoundaryAbsorption(fBoundaryAbsorptionCount);
248 run->IncHitsAboveThreshold(fPMTsAboveThreshold);
249
250 // If we have set the flag to save 'special' events, save here
251 if(fPhotonCount_Scint + fPhotonCount_Ceren < fDetector->GetSaveThreshold())
252 {
253 G4RunManager::GetRunManager()->rndmSaveThisEvent();
254 }
255}
256
257//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Definition of the LXeDetectorConstruction class.
Definition of the LXeEventAction class.
Definition of the LXeHistoManager class.
Definition of the LXePMTHit class.
G4THitsCollection< LXePMTHit > LXePMTHitsCollection
Definition LXePMTHit.hh:86
Definition of the LXeRun class.
Definition of the LXeScintHit class.
G4THitsCollection< LXeScintHit > LXeScintHitsCollection
Definition of the LXeTrajectory class.
~LXeEventAction() override
void EndOfEventAction(const G4Event *) override
void BeginOfEventAction(const G4Event *) override
G4ThreeVector fEWeightPos
G4ThreeVector fReconPos
LXeEventAction(const LXeDetectorConstruction *)
LXeEventMessenger * fEventMessenger
G4ThreeVector fPosMax
G4int fBoundaryAbsorptionCount
void IncHitCount(G4int count)
Definition LXeRun.hh:72

Applications | User Support | Publications | Collaboration