Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/optical/LXe/src/LXeEventAction.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

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