Geant4 Cross Reference |
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 B4/B4c/src/EventAction.cc 28 /// \brief Implementation of the B4c::EventAction class 29 30 #include "EventAction.hh" 31 32 #include "CalorHit.hh" 33 34 #include "G4AnalysisManager.hh" 35 #include "G4Event.hh" 36 #include "G4HCofThisEvent.hh" 37 #include "G4RunManager.hh" 38 #include "G4SDManager.hh" 39 #include "G4UnitsTable.hh" 40 41 #include <iomanip> 42 43 namespace B4c 44 { 45 46 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 47 48 CalorHitsCollection* EventAction::GetHitsCollection(G4int hcID, const G4Event* event) const 49 { 50 auto hitsCollection = static_cast<CalorHitsCollection*>(event->GetHCofThisEvent()->GetHC(hcID)); 51 52 if (!hitsCollection) { 53 G4ExceptionDescription msg; 54 msg << "Cannot access hitsCollection ID " << hcID; 55 G4Exception("EventAction::GetHitsCollection()", "MyCode0003", FatalException, msg); 56 } 57 58 return hitsCollection; 59 } 60 61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 62 63 void EventAction::PrintEventStatistics(G4double absoEdep, G4double absoTrackLength, 64 G4double gapEdep, G4double gapTrackLength) const 65 { 66 // print event statistics 67 G4cout << " Absorber: total energy: " << std::setw(7) << G4BestUnit(absoEdep, "Energy") 68 << " total track length: " << std::setw(7) << G4BestUnit(absoTrackLength, "Length") 69 << G4endl << " Gap: total energy: " << std::setw(7) << G4BestUnit(gapEdep, "Energy") 70 << " total track length: " << std::setw(7) << G4BestUnit(gapTrackLength, "Length") 71 << G4endl; 72 } 73 74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 75 76 void EventAction::BeginOfEventAction(const G4Event* /*event*/) {} 77 78 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 79 80 void EventAction::EndOfEventAction(const G4Event* event) 81 { 82 // Get hits collections IDs (only once) 83 if (fAbsHCID == -1) { 84 fAbsHCID = G4SDManager::GetSDMpointer()->GetCollectionID("AbsorberHitsCollection"); 85 fGapHCID = G4SDManager::GetSDMpointer()->GetCollectionID("GapHitsCollection"); 86 } 87 88 // Get hits collections 89 auto absoHC = GetHitsCollection(fAbsHCID, event); 90 auto gapHC = GetHitsCollection(fGapHCID, event); 91 92 // Get hit with total values 93 auto absoHit = (*absoHC)[absoHC->entries() - 1]; 94 auto gapHit = (*gapHC)[gapHC->entries() - 1]; 95 96 // Print per event (modulo n) 97 // 98 auto eventID = event->GetEventID(); 99 auto printModulo = G4RunManager::GetRunManager()->GetPrintProgress(); 100 if ((printModulo > 0) && (eventID % printModulo == 0)) { 101 PrintEventStatistics(absoHit->GetEdep(), absoHit->GetTrackLength(), gapHit->GetEdep(), 102 gapHit->GetTrackLength()); 103 G4cout << "--> End of event: " << eventID << "\n" << G4endl; 104 } 105 106 // Fill histograms, ntuple 107 // 108 109 // get analysis manager 110 auto analysisManager = G4AnalysisManager::Instance(); 111 112 // fill histograms 113 analysisManager->FillH1(0, absoHit->GetEdep()); 114 analysisManager->FillH1(1, gapHit->GetEdep()); 115 analysisManager->FillH1(2, absoHit->GetTrackLength()); 116 analysisManager->FillH1(3, gapHit->GetTrackLength()); 117 118 // fill ntuple 119 analysisManager->FillNtupleDColumn(0, absoHit->GetEdep()); 120 analysisManager->FillNtupleDColumn(1, gapHit->GetEdep()); 121 analysisManager->FillNtupleDColumn(2, absoHit->GetTrackLength()); 122 analysisManager->FillNtupleDColumn(3, gapHit->GetTrackLength()); 123 analysisManager->AddNtupleRow(); 124 } 125 126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 127 128 } // namespace B4c 129