Geant4  10.02.p03
B4dEventAction Class Reference

#include <B4dEventAction.hh>

Inheritance diagram for B4dEventAction:
Collaboration diagram for B4dEventAction:

Public Member Functions

 B4dEventAction ()
 
virtual ~B4dEventAction ()
 
virtual void BeginOfEventAction (const G4Event *event)
 
virtual void EndOfEventAction (const G4Event *event)
 
- Public Member Functions inherited from G4UserEventAction
 G4UserEventAction ()
 
virtual ~G4UserEventAction ()
 
virtual void SetEventManager (G4EventManager *value)
 

Private Member Functions

G4THitsMap< G4double > * GetHitsCollection (G4int hcID, const G4Event *event) const
 
G4double GetSum (G4THitsMap< G4double > *hitsMap) const
 
void PrintEventStatistics (G4double absoEdep, G4double absoTrackLength, G4double gapEdep, G4double gapTrackLength) const
 

Private Attributes

G4int fAbsoEdepHCID
 
G4int fGapEdepHCID
 
G4int fAbsoTrackLengthHCID
 
G4int fGapTrackLengthHCID
 

Additional Inherited Members

- Protected Attributes inherited from G4UserEventAction
G4EventManagerfpEventManager
 

Detailed Description

Event action class

In EndOfEventAction(), it prints the accumulated quantities of the energy deposit and track lengths of charged particles in Absober and Gap layers stored in the hits collections.

Definition at line 45 of file B4dEventAction.hh.

Constructor & Destructor Documentation

◆ B4dEventAction()

B4dEventAction::B4dEventAction ( )

Definition at line 45 of file B4dEventAction.cc.

47  fAbsoEdepHCID(-1),
48  fGapEdepHCID(-1),
51 {}
G4int fAbsoTrackLengthHCID
G4int fGapTrackLengthHCID

◆ ~B4dEventAction()

B4dEventAction::~B4dEventAction ( )
virtual

Definition at line 55 of file B4dEventAction.cc.

56 {}

Member Function Documentation

◆ BeginOfEventAction()

void B4dEventAction::BeginOfEventAction ( const G4Event event)
virtual

Reimplemented from G4UserEventAction.

Definition at line 113 of file B4dEventAction.cc.

114 {}

◆ EndOfEventAction()

void B4dEventAction::EndOfEventAction ( const G4Event event)
virtual

Reimplemented from G4UserEventAction.

Definition at line 118 of file B4dEventAction.cc.

119 {
120  // Get hist collections IDs
121  if ( fAbsoEdepHCID == -1 ) {
123  = G4SDManager::GetSDMpointer()->GetCollectionID("Absorber/Edep");
124  fGapEdepHCID
127  = G4SDManager::GetSDMpointer()->GetCollectionID("Absorber/TrackLength");
129  = G4SDManager::GetSDMpointer()->GetCollectionID("Gap/TrackLength");
130  }
131 
132  // Get sum values from hits collections
133  //
134  G4double absoEdep = GetSum(GetHitsCollection(fAbsoEdepHCID, event));
135  G4double gapEdep = GetSum(GetHitsCollection(fGapEdepHCID, event));
136 
137  G4double absoTrackLength
139  G4double gapTrackLength
141 
142  // get analysis manager
143  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
144 
145  // fill histograms
146  //
147  analysisManager->FillH1(1, absoEdep);
148  analysisManager->FillH1(2, gapEdep);
149  analysisManager->FillH1(3, absoTrackLength);
150  analysisManager->FillH1(4, gapTrackLength);
151 
152  // fill ntuple
153  //
154  analysisManager->FillNtupleDColumn(0, absoEdep);
155  analysisManager->FillNtupleDColumn(1, gapEdep);
156  analysisManager->FillNtupleDColumn(2, absoTrackLength);
157  analysisManager->FillNtupleDColumn(3, gapTrackLength);
158  analysisManager->AddNtupleRow();
159 
160  //print per event (modulo n)
161  //
162  G4int eventID = event->GetEventID();
164  if ( ( printModulo > 0 ) && ( eventID % printModulo == 0 ) ) {
165  G4cout << "---> End of event: " << eventID << G4endl;
166  PrintEventStatistics(absoEdep, absoTrackLength, gapEdep, gapTrackLength);
167  }
168 }
G4int GetPrintProgress()
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:135
G4THitsMap< G4double > * GetHitsCollection(G4int hcID, const G4Event *event) const
int G4int
Definition: G4Types.hh:78
void PrintEventStatistics(G4double absoEdep, G4double absoTrackLength, G4double gapEdep, G4double gapTrackLength) const
G4int fAbsoTrackLengthHCID
G4GLOB_DLL std::ostream G4cout
ExG4HbookAnalysisManager G4AnalysisManager
Definition: g4hbook_defs.hh:66
static G4RunManager * GetRunManager()
Definition: G4RunManager.cc:79
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
#define G4endl
Definition: G4ios.hh:61
G4int fGapTrackLengthHCID
G4double GetSum(G4THitsMap< G4double > *hitsMap) const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetHitsCollection()

G4THitsMap< G4double > * B4dEventAction::GetHitsCollection ( G4int  hcID,
const G4Event event 
) const
private

Definition at line 61 of file B4dEventAction.cc.

63 {
64  G4THitsMap<G4double>* hitsCollection
65  = static_cast<G4THitsMap<G4double>*>(
66  event->GetHCofThisEvent()->GetHC(hcID));
67 
68  if ( ! hitsCollection ) {
70  msg << "Cannot access hitsCollection ID " << hcID;
71  G4Exception("B4dEventAction::GetHitsCollection()",
72  "MyCode0003", FatalException, msg);
73  }
74 
75  return hitsCollection;
76 }
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetSum()

G4double B4dEventAction::GetSum ( G4THitsMap< G4double > *  hitsMap) const
private

Definition at line 80 of file B4dEventAction.cc.

81 {
82  G4double sumValue = 0;
83  std::map<G4int, G4double*>::iterator it;
84  for ( it = hitsMap->GetMap()->begin(); it != hitsMap->GetMap()->end(); it++) {
85  sumValue += *(it->second);
86  }
87  return sumValue;
88 }
std::map< G4int, T * > * GetMap() const
Definition: G4THitsMap.hh:68
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ PrintEventStatistics()

void B4dEventAction::PrintEventStatistics ( G4double  absoEdep,
G4double  absoTrackLength,
G4double  gapEdep,
G4double  gapTrackLength 
) const
private

Definition at line 92 of file B4dEventAction.cc.

95 {
96  // Print event statistics
97  //
98  G4cout
99  << " Absorber: total energy: "
100  << std::setw(7) << G4BestUnit(absoEdep, "Energy")
101  << " total track length: "
102  << std::setw(7) << G4BestUnit(absoTrackLength, "Length")
103  << G4endl
104  << " Gap: total energy: "
105  << std::setw(7) << G4BestUnit(gapEdep, "Energy")
106  << " total track length: "
107  << std::setw(7) << G4BestUnit(gapTrackLength, "Length")
108  << G4endl;
109 }
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
Here is the caller graph for this function:

Member Data Documentation

◆ fAbsoEdepHCID

G4int B4dEventAction::fAbsoEdepHCID
private

Definition at line 63 of file B4dEventAction.hh.

◆ fAbsoTrackLengthHCID

G4int B4dEventAction::fAbsoTrackLengthHCID
private

Definition at line 65 of file B4dEventAction.hh.

◆ fGapEdepHCID

G4int B4dEventAction::fGapEdepHCID
private

Definition at line 64 of file B4dEventAction.hh.

◆ fGapTrackLengthHCID

G4int B4dEventAction::fGapTrackLengthHCID
private

Definition at line 66 of file B4dEventAction.hh.


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