Geant4  10.02.p03
LXeEventAction Class Reference

#include <LXeEventAction.hh>

Inheritance diagram for LXeEventAction:
Collaboration diagram for LXeEventAction:

Public Member Functions

 LXeEventAction (LXeRecorderBase *)
 
virtual ~LXeEventAction ()
 
virtual void BeginOfEventAction (const G4Event *)
 
virtual void EndOfEventAction (const G4Event *)
 
void SetSaveThreshold (G4int)
 
void SetEventVerbose (G4int v)
 
void SetPMTThreshold (G4int t)
 
void SetForceDrawPhotons (G4bool b)
 
void SetForceDrawNoPhotons (G4bool b)
 
- Public Member Functions inherited from G4UserEventAction
 G4UserEventAction ()
 
virtual ~G4UserEventAction ()
 
virtual void SetEventManager (G4EventManager *value)
 

Private Attributes

LXeRecorderBasefRecorder
 
LXeEventMessengerfEventMessenger
 
G4int fSaveThreshold
 
G4int fScintCollID
 
G4int fPMTCollID
 
G4int fVerbose
 
G4int fPMTThreshold
 
G4bool fForcedrawphotons
 
G4bool fForcenophotons
 

Additional Inherited Members

- Protected Attributes inherited from G4UserEventAction
G4EventManagerfpEventManager
 

Detailed Description

Definition at line 43 of file LXeEventAction.hh.

Constructor & Destructor Documentation

◆ LXeEventAction()

LXeEventAction::LXeEventAction ( LXeRecorderBase r)

Definition at line 53 of file LXeEventAction.cc.

56 {
58 }
LXeRecorderBase * fRecorder
LXeEventMessenger * fEventMessenger
G4bool fForcedrawphotons
G4bool fForcenophotons

◆ ~LXeEventAction()

LXeEventAction::~LXeEventAction ( )
virtual

Definition at line 62 of file LXeEventAction.cc.

62 {}

Member Function Documentation

◆ BeginOfEventAction()

void LXeEventAction::BeginOfEventAction ( const G4Event anEvent)
virtual

Reimplemented from G4UserEventAction.

Definition at line 66 of file LXeEventAction.cc.

66  {
67 
68  //New event, add the user information object
71 
73  if(fScintCollID<0)
74  fScintCollID=SDman->GetCollectionID("scintCollection");
75  if(fPMTCollID<0)
76  fPMTCollID=SDman->GetCollectionID("pmtHitCollection");
77 
79 }
LXeRecorderBase * fRecorder
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:135
void SetUserInformation(G4VUserEventInformation *anInfo)
virtual void RecordBeginOfEvent(const G4Event *)
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
static G4EventManager * GetEventManager()
Here is the call graph for this function:

◆ EndOfEventAction()

void LXeEventAction::EndOfEventAction ( const G4Event anEvent)
virtual

Reimplemented from G4UserEventAction.

Definition at line 83 of file LXeEventAction.cc.

83  {
84 
85  LXeUserEventInformation* eventInformation
87 
88  G4TrajectoryContainer* trajectoryContainer=anEvent->GetTrajectoryContainer();
89 
90  G4int n_trajectories = 0;
91  if (trajectoryContainer) n_trajectories = trajectoryContainer->entries();
92 
93  // extract the trajectories and draw them
95  for (G4int i=0; i<n_trajectories; i++){
97  ((*(anEvent->GetTrajectoryContainer()))[i]);
98  if(trj->GetParticleName()=="opticalphoton"){
101  }
102  trj->DrawTrajectory();
103  }
104  }
105 
106  LXeScintHitsCollection* scintHC = 0;
107  LXePMTHitsCollection* pmtHC = 0;
108  G4HCofThisEvent* hitsCE = anEvent->GetHCofThisEvent();
109 
110  //Get the hit collections
111  if(hitsCE){
112  if(fScintCollID>=0)scintHC = (LXeScintHitsCollection*)(hitsCE->GetHC(fScintCollID));
113  if(fPMTCollID>=0)pmtHC = (LXePMTHitsCollection*)(hitsCE->GetHC(fPMTCollID));
114  }
115 
116  //Hits in scintillator
117  if(scintHC){
118  int n_hit = scintHC->entries();
119  G4ThreeVector eWeightPos(0.);
120  G4double edep;
121  G4double edepMax=0;
122 
123  for(int i=0;i<n_hit;i++){ //gather info on hits in scintillator
124  edep=(*scintHC)[i]->GetEdep();
125  eventInformation->IncEDep(edep); //sum up the edep
126  eWeightPos += (*scintHC)[i]->GetPos()*edep;//calculate energy weighted pos
127  if(edep>edepMax){
128  edepMax=edep;//store max energy deposit
129  G4ThreeVector posMax=(*scintHC)[i]->GetPos();
130  eventInformation->SetPosMax(posMax,edep);
131  }
132  }
133  if(eventInformation->GetEDep()==0.){
134  if(fVerbose>0)G4cout<<"No hits in the scintillator this event."<<G4endl;
135  }
136  else{
137  //Finish calculation of energy weighted position
138  eWeightPos/=eventInformation->GetEDep();
139  eventInformation->SetEWeightPos(eWeightPos);
140  if(fVerbose>0){
141  G4cout << "\tEnergy weighted position of hits in LXe : "
142  << eWeightPos/mm << G4endl;
143  }
144  }
145  if(fVerbose>0){
146  G4cout << "\tTotal energy deposition in scintillator : "
147  << eventInformation->GetEDep() / keV << " (keV)" << G4endl;
148  }
149  }
150 
151  if(pmtHC){
152  G4ThreeVector reconPos(0.,0.,0.);
153  G4int pmts=pmtHC->entries();
154  //Gather info from all PMTs
155  for(G4int i=0;i<pmts;i++){
156  eventInformation->IncHitCount((*pmtHC)[i]->GetPhotonCount());
157  reconPos+=(*pmtHC)[i]->GetPMTPos()*(*pmtHC)[i]->GetPhotonCount();
158  if((*pmtHC)[i]->GetPhotonCount()>=fPMTThreshold){
159  eventInformation->IncPMTSAboveThreshold();
160  }
161  else{//wasnt above the threshold, turn it back off
162  (*pmtHC)[i]->SetDrawit(false);
163  }
164  }
165 
166  if(eventInformation->GetHitCount()>0){//dont bother unless there were hits
167  reconPos/=eventInformation->GetHitCount();
168  if(fVerbose>0){
169  G4cout << "\tReconstructed position of hits in LXe : "
170  << reconPos/mm << G4endl;
171  }
172  eventInformation->SetReconPos(reconPos);
173  }
174  pmtHC->DrawAllHits();
175  }
176 
177  if(fVerbose>0){
178  //End of event output. later to be controlled by a verbose level
179  G4cout << "\tNumber of photons that hit PMTs in this event : "
180  << eventInformation->GetHitCount() << G4endl;
181  G4cout << "\tNumber of PMTs above threshold("<<fPMTThreshold<<") : "
182  << eventInformation->GetPMTSAboveThreshold() << G4endl;
183  G4cout << "\tNumber of photons produced by scintillation in this event : "
184  << eventInformation->GetPhotonCount_Scint() << G4endl;
185  G4cout << "\tNumber of photons produced by cerenkov in this event : "
186  << eventInformation->GetPhotonCount_Ceren() << G4endl;
187  G4cout << "\tNumber of photons absorbed (OpAbsorption) in this event : "
188  << eventInformation->GetAbsorptionCount() << G4endl;
189  G4cout << "\tNumber of photons absorbed at boundaries (OpBoundary) in "
190  << "this event : " << eventInformation->GetBoundaryAbsorptionCount()
191  << G4endl;
192  G4cout << "Unacounted for photons in this event : "
193  << (eventInformation->GetPhotonCount_Scint() +
194  eventInformation->GetPhotonCount_Ceren() -
195  eventInformation->GetAbsorptionCount() -
196  eventInformation->GetHitCount() -
197  eventInformation->GetBoundaryAbsorptionCount())
198  << G4endl;
199  }
200  //If we have set the flag to save 'special' events, save here
201  if(fSaveThreshold&&eventInformation->GetPhotonCount() <= fSaveThreshold)
203 
204  if(fRecorder)fRecorder->RecordEndOfEvent(anEvent);
205 }
LXeRecorderBase * fRecorder
static G4VVisManager * GetConcreteInstance()
int G4int
Definition: G4Types.hh:78
Double_t edep
G4VUserEventInformation * GetUserInformation() const
Definition: G4Event.hh:199
G4GLOB_DLL std::ostream G4cout
virtual void DrawTrajectory() const
virtual void rndmSaveThisEvent()
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
G4TrajectoryContainer * GetTrajectoryContainer() const
Definition: G4Event.hh:189
void SetReconPos(const G4ThreeVector &p)
G4String GetParticleName() const
Definition: G4Trajectory.hh:99
G4HCofThisEvent * GetHCofThisEvent() const
Definition: G4Event.hh:185
#define G4endl
Definition: G4ios.hh:61
G4bool fForcedrawphotons
G4bool fForcenophotons
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
void SetForceNoDrawTrajectory(G4bool b)
static const double mm
Definition: G4SIunits.hh:114
virtual void RecordEndOfEvent(const G4Event *)
void SetEWeightPos(const G4ThreeVector &p)
void SetForceDrawTrajectory(G4bool b)
void SetPosMax(const G4ThreeVector &p, G4double edep)
Here is the call graph for this function:

◆ SetEventVerbose()

void LXeEventAction::SetEventVerbose ( G4int  v)
inline

Definition at line 57 of file LXeEventAction.hh.

Here is the caller graph for this function:

◆ SetForceDrawNoPhotons()

void LXeEventAction::SetForceDrawNoPhotons ( G4bool  b)
inline

Definition at line 62 of file LXeEventAction.hh.

Here is the caller graph for this function:

◆ SetForceDrawPhotons()

void LXeEventAction::SetForceDrawPhotons ( G4bool  b)
inline

Definition at line 61 of file LXeEventAction.hh.

Here is the caller graph for this function:

◆ SetPMTThreshold()

void LXeEventAction::SetPMTThreshold ( G4int  t)
inline

Definition at line 59 of file LXeEventAction.hh.

59 {fPMTThreshold=t;}
Here is the caller graph for this function:

◆ SetSaveThreshold()

void LXeEventAction::SetSaveThreshold ( G4int  save)

Definition at line 209 of file LXeEventAction.cc.

209  {
210 /*Sets the save threshold for the random number seed. If the number of photons
211 generated in an event is lower than this, then save the seed for this event
212 in a file called run###evt###.rndm
213 */
214  fSaveThreshold=save;
217  // G4UImanager::GetUIpointer()->ApplyCommand("/random/setSavingFlag 1");
218 }
void SetRandomNumberStore(G4bool flag)
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
void SetRandomNumberStoreDir(const G4String &dir)
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ fEventMessenger

LXeEventMessenger* LXeEventAction::fEventMessenger
private

Definition at line 67 of file LXeEventAction.hh.

◆ fForcedrawphotons

G4bool LXeEventAction::fForcedrawphotons
private

Definition at line 78 of file LXeEventAction.hh.

◆ fForcenophotons

G4bool LXeEventAction::fForcenophotons
private

Definition at line 79 of file LXeEventAction.hh.

◆ fPMTCollID

G4int LXeEventAction::fPMTCollID
private

Definition at line 72 of file LXeEventAction.hh.

◆ fPMTThreshold

G4int LXeEventAction::fPMTThreshold
private

Definition at line 76 of file LXeEventAction.hh.

◆ fRecorder

LXeRecorderBase* LXeEventAction::fRecorder
private

Definition at line 66 of file LXeEventAction.hh.

◆ fSaveThreshold

G4int LXeEventAction::fSaveThreshold
private

Definition at line 69 of file LXeEventAction.hh.

◆ fScintCollID

G4int LXeEventAction::fScintCollID
private

Definition at line 71 of file LXeEventAction.hh.

◆ fVerbose

G4int LXeEventAction::fVerbose
private

Definition at line 74 of file LXeEventAction.hh.


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