Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
LXeEventAction Class Reference

#include <Doxymodules_optical.h>

Inheritance diagram for LXeEventAction:
G4UserEventAction

Public Member Functions

 LXeEventAction (const LXeDetectorConstruction *)
 
 ~LXeEventAction () override
 
void BeginOfEventAction (const G4Event *) override
 
void EndOfEventAction (const G4Event *) override
 
void SetEventVerbose (G4int v)
 
void SetPMTThreshold (G4int t)
 
void SetForceDrawPhotons (G4bool b)
 
void SetForceDrawNoPhotons (G4bool b)
 
void IncPhotonCount_Scint ()
 
void IncPhotonCount_Ceren ()
 
void IncEDep (G4double dep)
 
void IncAbsorption ()
 
void IncBoundaryAbsorption ()
 
void IncHitCount (G4int i=1)
 
void SetEWeightPos (const G4ThreeVector &p)
 
void SetReconPos (const G4ThreeVector &p)
 
void SetConvPos (const G4ThreeVector &p)
 
void SetPosMax (const G4ThreeVector &p, G4double edep)
 
G4int GetPhotonCount_Scint () const
 
G4int GetPhotonCount_Ceren () const
 
G4int GetHitCount () const
 
G4double GetEDep () const
 
G4int GetAbsorptionCount () const
 
G4int GetBoundaryAbsorptionCount () const
 
G4ThreeVector GetEWeightPos ()
 
G4ThreeVector GetReconPos ()
 
G4ThreeVector GetConvPos ()
 
G4ThreeVector GetPosMax ()
 
G4double GetEDepMax ()
 
G4double IsConvPosSet ()
 
G4int GetPhotonCount ()
 
void IncPMTSAboveThreshold ()
 
G4int GetPMTSAboveThreshold ()
 

Private Attributes

LXeEventMessengerfEventMessenger = nullptr
 
const LXeDetectorConstructionfDetector = nullptr
 
G4int fScintCollID = -1
 
G4int fPMTCollID = -1
 
G4int fVerbose = 0
 
G4int fPMTThreshold = 1
 
G4bool fForcedrawphotons = false
 
G4bool fForcenophotons = false
 
G4int fHitCount = 0
 
G4int fPhotonCount_Scint = 0
 
G4int fPhotonCount_Ceren = 0
 
G4int fAbsorptionCount = 0
 
G4int fBoundaryAbsorptionCount = 0
 
G4double fTotE = 0.
 
G4ThreeVector fEWeightPos
 
G4ThreeVector fReconPos
 
G4ThreeVector fConvPos
 
G4bool fConvPosSet = false
 
G4ThreeVector fPosMax
 
G4double fEdepMax = 0.
 
G4int fPMTsAboveThreshold = 0
 

Detailed Description

Definition at line 51 of file Doxymodules_optical.h.

Constructor & Destructor Documentation

◆ LXeEventAction()

LXeEventAction::LXeEventAction ( const LXeDetectorConstruction det)

Definition at line 53 of file LXeEventAction.cc.

54 : fDetector(det)
55{
57}
const LXeDetectorConstruction * fDetector
LXeEventMessenger * fEventMessenger

◆ ~LXeEventAction()

LXeEventAction::~LXeEventAction ( )
override

Definition at line 61 of file LXeEventAction.cc.

61{ delete fEventMessenger; }

Member Function Documentation

◆ BeginOfEventAction()

void LXeEventAction::BeginOfEventAction ( const G4Event )
override

Definition at line 65 of file LXeEventAction.cc.

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}
G4int fBoundaryAbsorptionCount

◆ EndOfEventAction()

void LXeEventAction::EndOfEventAction ( const G4Event anEvent)
override

Definition at line 88 of file LXeEventAction.cc.

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}
G4THitsCollection< LXePMTHit > LXePMTHitsCollection
Definition LXePMTHit.hh:86
G4THitsCollection< LXeScintHit > LXeScintHitsCollection
G4ThreeVector fEWeightPos
G4ThreeVector fReconPos
G4ThreeVector fPosMax
void IncHitCount(G4int count)
Definition LXeRun.hh:72

◆ SetEventVerbose()

void LXeEventAction::SetEventVerbose ( G4int  v)
inline

Definition at line 52 of file LXeEventAction.hh.

52{ fVerbose = v; }

◆ SetPMTThreshold()

void LXeEventAction::SetPMTThreshold ( G4int  t)
inline

Definition at line 54 of file LXeEventAction.hh.

54{ fPMTThreshold = t; }

◆ SetForceDrawPhotons()

void LXeEventAction::SetForceDrawPhotons ( G4bool  b)
inline

Definition at line 56 of file LXeEventAction.hh.

56{ fForcedrawphotons = b; }

◆ SetForceDrawNoPhotons()

void LXeEventAction::SetForceDrawNoPhotons ( G4bool  b)
inline

Definition at line 57 of file LXeEventAction.hh.

57{ fForcenophotons = b; }

◆ IncPhotonCount_Scint()

void LXeEventAction::IncPhotonCount_Scint ( )
inline

Definition at line 59 of file LXeEventAction.hh.

◆ IncPhotonCount_Ceren()

void LXeEventAction::IncPhotonCount_Ceren ( )
inline

Definition at line 60 of file LXeEventAction.hh.

◆ IncEDep()

void LXeEventAction::IncEDep ( G4double  dep)
inline

Definition at line 61 of file LXeEventAction.hh.

61{ fTotE += dep; }

◆ IncAbsorption()

void LXeEventAction::IncAbsorption ( )
inline

Definition at line 62 of file LXeEventAction.hh.

◆ IncBoundaryAbsorption()

void LXeEventAction::IncBoundaryAbsorption ( )
inline

Definition at line 63 of file LXeEventAction.hh.

◆ IncHitCount()

void LXeEventAction::IncHitCount ( G4int  i = 1)
inline

Definition at line 64 of file LXeEventAction.hh.

64{ fHitCount += i; }

◆ SetEWeightPos()

void LXeEventAction::SetEWeightPos ( const G4ThreeVector &  p)
inline

Definition at line 66 of file LXeEventAction.hh.

66{ fEWeightPos = p; }

◆ SetReconPos()

void LXeEventAction::SetReconPos ( const G4ThreeVector &  p)
inline

Definition at line 67 of file LXeEventAction.hh.

67{ fReconPos = p; }

◆ SetConvPos()

void LXeEventAction::SetConvPos ( const G4ThreeVector &  p)
inline

Definition at line 68 of file LXeEventAction.hh.

69 {
70 fConvPos = p;
71 fConvPosSet = true;
72 }
G4ThreeVector fConvPos

◆ SetPosMax()

void LXeEventAction::SetPosMax ( const G4ThreeVector &  p,
G4double  edep 
)
inline

Definition at line 73 of file LXeEventAction.hh.

74 {
75 fPosMax = p;
76 fEdepMax = edep;
77 }

◆ GetPhotonCount_Scint()

G4int LXeEventAction::GetPhotonCount_Scint ( ) const
inline

Definition at line 79 of file LXeEventAction.hh.

79{ return fPhotonCount_Scint; }

◆ GetPhotonCount_Ceren()

G4int LXeEventAction::GetPhotonCount_Ceren ( ) const
inline

Definition at line 80 of file LXeEventAction.hh.

80{ return fPhotonCount_Ceren; }

◆ GetHitCount()

G4int LXeEventAction::GetHitCount ( ) const
inline

Definition at line 81 of file LXeEventAction.hh.

81{ return fHitCount; }

◆ GetEDep()

G4double LXeEventAction::GetEDep ( ) const
inline

Definition at line 82 of file LXeEventAction.hh.

82{ return fTotE; }

◆ GetAbsorptionCount()

G4int LXeEventAction::GetAbsorptionCount ( ) const
inline

Definition at line 83 of file LXeEventAction.hh.

83{ return fAbsorptionCount; }

◆ GetBoundaryAbsorptionCount()

G4int LXeEventAction::GetBoundaryAbsorptionCount ( ) const
inline

Definition at line 84 of file LXeEventAction.hh.

◆ GetEWeightPos()

G4ThreeVector LXeEventAction::GetEWeightPos ( )
inline

Definition at line 86 of file LXeEventAction.hh.

86{ return fEWeightPos; }

◆ GetReconPos()

G4ThreeVector LXeEventAction::GetReconPos ( )
inline

Definition at line 87 of file LXeEventAction.hh.

87{ return fReconPos; }

◆ GetConvPos()

G4ThreeVector LXeEventAction::GetConvPos ( )
inline

Definition at line 88 of file LXeEventAction.hh.

88{ return fConvPos; }

◆ GetPosMax()

G4ThreeVector LXeEventAction::GetPosMax ( )
inline

Definition at line 89 of file LXeEventAction.hh.

89{ return fPosMax; }

◆ GetEDepMax()

G4double LXeEventAction::GetEDepMax ( )
inline

Definition at line 90 of file LXeEventAction.hh.

90{ return fEdepMax; }

◆ IsConvPosSet()

G4double LXeEventAction::IsConvPosSet ( )
inline

Definition at line 91 of file LXeEventAction.hh.

91{ return fConvPosSet; }

◆ GetPhotonCount()

G4int LXeEventAction::GetPhotonCount ( )
inline

Definition at line 94 of file LXeEventAction.hh.

◆ IncPMTSAboveThreshold()

void LXeEventAction::IncPMTSAboveThreshold ( )
inline

Definition at line 96 of file LXeEventAction.hh.

◆ GetPMTSAboveThreshold()

G4int LXeEventAction::GetPMTSAboveThreshold ( )
inline

Definition at line 97 of file LXeEventAction.hh.

97{ return fPMTsAboveThreshold; }

Member Data Documentation

◆ fEventMessenger

LXeEventMessenger* LXeEventAction::fEventMessenger = nullptr
private

Definition at line 100 of file LXeEventAction.hh.

◆ fDetector

const LXeDetectorConstruction* LXeEventAction::fDetector = nullptr
private

Definition at line 101 of file LXeEventAction.hh.

◆ fScintCollID

G4int LXeEventAction::fScintCollID = -1
private

Definition at line 103 of file LXeEventAction.hh.

◆ fPMTCollID

G4int LXeEventAction::fPMTCollID = -1
private

Definition at line 104 of file LXeEventAction.hh.

◆ fVerbose

G4int LXeEventAction::fVerbose = 0
private

Definition at line 106 of file LXeEventAction.hh.

◆ fPMTThreshold

G4int LXeEventAction::fPMTThreshold = 1
private

Definition at line 108 of file LXeEventAction.hh.

◆ fForcedrawphotons

G4bool LXeEventAction::fForcedrawphotons = false
private

Definition at line 110 of file LXeEventAction.hh.

◆ fForcenophotons

G4bool LXeEventAction::fForcenophotons = false
private

Definition at line 111 of file LXeEventAction.hh.

◆ fHitCount

G4int LXeEventAction::fHitCount = 0
private

Definition at line 113 of file LXeEventAction.hh.

◆ fPhotonCount_Scint

G4int LXeEventAction::fPhotonCount_Scint = 0
private

Definition at line 114 of file LXeEventAction.hh.

◆ fPhotonCount_Ceren

G4int LXeEventAction::fPhotonCount_Ceren = 0
private

Definition at line 115 of file LXeEventAction.hh.

◆ fAbsorptionCount

G4int LXeEventAction::fAbsorptionCount = 0
private

Definition at line 116 of file LXeEventAction.hh.

◆ fBoundaryAbsorptionCount

G4int LXeEventAction::fBoundaryAbsorptionCount = 0
private

Definition at line 117 of file LXeEventAction.hh.

◆ fTotE

G4double LXeEventAction::fTotE = 0.
private

Definition at line 119 of file LXeEventAction.hh.

◆ fEWeightPos

G4ThreeVector LXeEventAction::fEWeightPos
private

Definition at line 123 of file LXeEventAction.hh.

◆ fReconPos

G4ThreeVector LXeEventAction::fReconPos
private

Definition at line 124 of file LXeEventAction.hh.

◆ fConvPos

G4ThreeVector LXeEventAction::fConvPos
private

Definition at line 125 of file LXeEventAction.hh.

◆ fConvPosSet

G4bool LXeEventAction::fConvPosSet = false
private

Definition at line 126 of file LXeEventAction.hh.

◆ fPosMax

G4ThreeVector LXeEventAction::fPosMax
private

Definition at line 127 of file LXeEventAction.hh.

◆ fEdepMax

G4double LXeEventAction::fEdepMax = 0.
private

Definition at line 128 of file LXeEventAction.hh.

◆ fPMTsAboveThreshold

G4int LXeEventAction::fPMTsAboveThreshold = 0
private

Definition at line 130 of file LXeEventAction.hh.


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

Applications | User Support | Publications | Collaboration