Geant4 examples
basic/B4/B4d/src/B4dEventAction.cc
Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00030 
00031 #include "B4dEventAction.hh"
00032 #include "B4dEventActionMessenger.hh"
00033 #include "B4Analysis.hh"
00034 
00035 #include "G4RunManager.hh"
00036 #include "G4Event.hh"
00037 #include "G4SDManager.hh"
00038 #include "G4HCofThisEvent.hh"
00039 #include "G4UnitsTable.hh"
00040 
00041 #include "Randomize.hh"
00042 #include <iomanip>
00043 
00044 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00045 
00046 B4dEventAction::B4dEventAction()
00047  : G4UserEventAction(),
00048    fMessenger(this),
00049    fPrintModulo(1)
00050 {
00051 }
00052 
00053 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00054 
00055 B4dEventAction::~B4dEventAction()
00056 {
00057 }
00058 
00059 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00060 
00061 G4THitsMap<G4double>* 
00062 B4dEventAction::GetHitsCollection(const G4String& hcName,
00063                                   const G4Event* event) const
00064 {
00065   G4int hcID 
00066     = G4SDManager::GetSDMpointer()->GetCollectionID(hcName);
00067   G4THitsMap<G4double>* hitsCollection 
00068     = static_cast<G4THitsMap<G4double>*>(
00069         event->GetHCofThisEvent()->GetHC(hcID));
00070   
00071   if ( ! hitsCollection ) {
00072     G4cerr << "Cannot access hitsCollection " << hcName << G4endl;
00073     exit(1);
00074   }         
00075 
00076   return hitsCollection;
00077 }    
00078 
00079 G4double B4dEventAction::GetSum(G4THitsMap<G4double>* hitsMap) const
00080 {
00081   G4double sumValue = 0;
00082   std::map<G4int, G4double*>::iterator it;
00083   for ( it = hitsMap->GetMap()->begin(); it != hitsMap->GetMap()->end(); it++) {
00084     sumValue += *(it->second);
00085   }
00086   return sumValue;  
00087 }  
00088 
00089 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00090 
00091 void B4dEventAction::PrintEventStatistics(
00092                             G4double absoEdep, G4double absoTrackLength,
00093                             G4double gapEdep, G4double gapTrackLength) const
00094 {
00095   // Print event statistics
00096   //
00097   G4cout
00098      << "   Absorber: total energy: " 
00099      << std::setw(7) << G4BestUnit(absoEdep, "Energy")
00100      << "       total track length: " 
00101      << std::setw(7) << G4BestUnit(absoTrackLength, "Length")
00102      << G4endl
00103      << "        Gap: total energy: " 
00104      << std::setw(7) << G4BestUnit(gapEdep, "Energy")
00105      << "       total track length: " 
00106      << std::setw(7) << G4BestUnit(gapTrackLength, "Length")
00107      << G4endl;
00108 }
00109 
00110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00111 
00112 void B4dEventAction::BeginOfEventAction(const G4Event* event)
00113 {  
00114 
00115   G4int eventID = event->GetEventID();
00116   if ( eventID % fPrintModulo == 0 )  { 
00117     G4cout << "\n---> Begin of event: " << eventID << G4endl;
00118     //CLHEP::HepRandom::showEngineStatus();
00119   }
00120 }
00121 
00122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00123 
00124 void B4dEventAction::EndOfEventAction(const G4Event* event)
00125 {  
00126   // Get sum value from hits collections
00127   //
00128   G4double absoEdep 
00129     = GetSum(GetHitsCollection("Absorber/Edep", event));
00130 
00131   G4double gapEdep 
00132     = GetSum(GetHitsCollection("Gap/Edep", event));
00133 
00134   G4double absoTrackLength 
00135     = GetSum(GetHitsCollection("Absorber/TrackLength", event));
00136 
00137   G4double gapTrackLength 
00138     = GetSum(GetHitsCollection("Gap/TrackLength", event));
00139 
00140   // get analysis manager
00141   G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
00142 
00143   // fill histograms
00144   //  
00145   analysisManager->FillH1(1, absoEdep);
00146   analysisManager->FillH1(2, gapEdep);
00147   analysisManager->FillH1(3, absoTrackLength);
00148   analysisManager->FillH1(4, gapTrackLength);
00149   
00150   // fill ntuple
00151   //
00152   analysisManager->FillNtupleDColumn(0, absoEdep);
00153   analysisManager->FillNtupleDColumn(1, gapEdep);
00154   analysisManager->FillNtupleDColumn(2, absoTrackLength);
00155   analysisManager->FillNtupleDColumn(3, gapTrackLength);
00156   analysisManager->AddNtupleRow();  
00157   
00158   //print per event (modulo n)
00159   //
00160   G4int eventID = event->GetEventID();
00161   if ( eventID % fPrintModulo == 0) {
00162     G4cout << "---> End of event: " << eventID << G4endl;     
00163 
00164     PrintEventStatistics(absoEdep, absoTrackLength, gapEdep, gapTrackLength);
00165   }
00166 }  
00167 
00168 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines