Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ExGflashEventAction Class Reference

#include <ExGflashEventAction.hh>

Inheritance diagram for ExGflashEventAction:
Collaboration diagram for ExGflashEventAction:

Public Member Functions

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

Additional Inherited Members

- Protected Attributes inherited from G4UserEventAction
G4EventManagerfpEventManager
 

Detailed Description

Definition at line 38 of file ExGflashEventAction.hh.

Constructor & Destructor Documentation

ExGflashEventAction::ExGflashEventAction ( )

Definition at line 51 of file ExGflashEventAction.cc.

53  fNevent(0),fDtime(0.0),fCalorimeterCollectionId(-1)
54 {
55 }
ExGflashEventAction::~ExGflashEventAction ( )

Definition at line 59 of file ExGflashEventAction.cc.

60 {
61  if ( fNevent > 0 ) {
62  G4cout << "Internal Real Elapsed Time /event is: "<< fDtime /fNevent<< G4endl;
63  }
64 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Member Function Documentation

void ExGflashEventAction::BeginOfEventAction ( const G4Event evt)
virtual

Reimplemented from G4UserEventAction.

Definition at line 68 of file ExGflashEventAction.cc.

69 {
70  fTimerIntern.Start();
71  G4cout<<" ------ Start ExGflashEventAction ----- "<<G4endl;
72  fNevent=evt->GetEventID();
73  G4cout<<" Start generating event Nr "<<fNevent<<G4endl<<G4endl;
74 }
G4int GetEventID() const
Definition: G4Event.hh:151
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void Start()

Here is the call graph for this function:

void ExGflashEventAction::EndOfEventAction ( const G4Event evt)
virtual

Hits in sensitive Detector

Computing (x,y,z) of vertex of initial particles

For all Hits in sensitive detector

Reimplemented from G4UserEventAction.

Definition at line 78 of file ExGflashEventAction.cc.

79 {
80  fTimerIntern.Stop();
81  G4cout << G4endl;
82  G4cout << "******************************************";
83  G4cout << G4endl;
84  G4cout << "Internal Real Elapsed Time is: "<< fTimerIntern.GetRealElapsed();
85  G4cout << G4endl;
86  G4cout << "Internal System Elapsed Time: " << fTimerIntern.GetSystemElapsed();
87  G4cout << G4endl;
88  G4cout << "Internal GetUserElapsed Time: " << fTimerIntern.GetUserElapsed();
89  G4cout << G4endl;
90  G4cout << "******************************************"<< G4endl;
91  fDtime+=fTimerIntern.GetRealElapsed();
92  G4cout<<" ------ ExGflashEventAction::End of event nr. "<<fNevent<<" -----"<< G4endl;
93 
95  G4String colNam;
96  fCalorimeterCollectionId=SDman->GetCollectionID(colNam="ExGflashCollection");
97  if (fCalorimeterCollectionId<0) return;
98  G4HCofThisEvent * HCE = evt->GetHCofThisEvent();
99  ExGflashHitsCollection* THC = 0;
100  G4double totE = 0;
101  // Read out of the crysta ECAL
102  THC=(ExGflashHitsCollection *)(HCE->GetHC(fCalorimeterCollectionId));
103  if (THC)
104  {
106  int n_hit = THC->entries();
107  G4cout<<" " << n_hit<< " hits are stored in ExGflashHitsCollection "<<G4endl;
108  G4PrimaryVertex* pvertex=evt->GetPrimaryVertex();
110  G4ThreeVector vtx=pvertex->GetPosition();
111  G4PrimaryParticle* pparticle=pvertex->GetPrimary();
112  // direction of the Shower
113  G4ThreeVector mom=pparticle->GetMomentum()/pparticle->GetMomentum().mag();
114 
115  //@@@ ExGflashEventAction: Magicnumber
116  G4double energyincrystal[100];
117  for (int i=0;i<100;i++) energyincrystal[i]=0.;
118 
119  //@@@ ExGflashEventAction: Magicnumber
121  for (int i=0;i<n_hit;i++)
122  {
123  G4double estep = (*THC)[i]->GetEdep()/GeV;
124  if (estep >0.0)
125  {
126  totE += (*THC)[i]->GetEdep()/GeV;
127  G4int num=(*THC)[i]->GetCrystalNum();
128 
129  energyincrystal[num]+=(*THC)[i]->GetEdep()/GeV;
130  //G4cout << num << G4endl;
131  // G4cout << " Crystal Nummer " << (*THC)[i]->GetCrystalNum() << G4endl;
132  // G4cout << (*THC)[i]->GetCrystalNum() /10 <<
133  // " "<<(*THC)[i]->GetCrystalNum()%10 << G4endl;
134 
135  G4ThreeVector hitpos=(*THC)[i]->GetPos();
136  G4ThreeVector l (hitpos.x(), hitpos.y(), hitpos.z());
137  // distance from shower start
138  l = l - vtx;
139  // projection on shower axis = longitudinal profile
140  G4ThreeVector longitudinal = l;
141  // shower profiles (Radial)
142  G4ThreeVector radial = vtx.cross(l);
143  }
144  }
145  G4double max = 0;
146  G4int index = 0;
147  //Find crystal with maximum energy
148  for (int i=0;i<100;i++)
149  {
150  //G4cout << i <<" " << energyincrystal[i] << G4endl;
151  if (max <energyincrystal[i])
152  {
153  max=energyincrystal[i];
154  index = i;
155  }
156  }
157  //G4cout << index <<" NMAX " << index << G4endl;
158 
159  // 3x3 matrix of crystals around the crystal with the maximum energy deposit
160  G4double e3x3 =
161  energyincrystal[index]+energyincrystal[index+1]+energyincrystal[index-1]+
162  energyincrystal[index-10]+energyincrystal[index-9]+energyincrystal[index-11]+
163  energyincrystal[index+10]+energyincrystal[index+11]+energyincrystal[index+9];
164 
165  // 5x5 matrix of crystals around the crystal with the maximum energy deposit
166  G4double e5x5 =
167  energyincrystal[index]+energyincrystal[index+1]+energyincrystal[index-1]+
168  energyincrystal[index+2]+energyincrystal[index-2]+
169  energyincrystal[index-10]+energyincrystal[index-9]+energyincrystal[index-11]+
170  energyincrystal[index-8]+energyincrystal[index-12]+
171  energyincrystal[index+10]+energyincrystal[index+11]+energyincrystal[index+9]+
172  energyincrystal[index+12]+energyincrystal[index+8];
173 
174  G4cout << " e1 " << energyincrystal[index]
175  << " e3x3 " << e3x3<< " GeV e5x5 " <<e5x5 <<G4endl;
176  }
177 
178  G4cout << " Total energy deposited in the calorimeter: " << totE << " (GeV)" << G4endl;
179  G4TrajectoryContainer * trajectoryContainer = evt->GetTrajectoryContainer();
180  G4int n_trajectories = 0;
181  if(trajectoryContainer){ n_trajectories = trajectoryContainer->entries(); }
182  G4cout << " " << n_trajectories << " trajectories stored in this event." << G4endl;
183 
184 }
G4ThreeVector GetMomentum() const
G4double GetSystemElapsed() const
Definition: G4Timer.cc:119
G4VHitsCollection * GetHC(G4int i)
G4ThreeVector GetPosition() const
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:135
double x() const
G4int entries() const
int G4int
Definition: G4Types.hh:78
G4TrajectoryContainer * GetTrajectoryContainer() const
Definition: G4Event.hh:189
double z() const
G4PrimaryParticle * GetPrimary(G4int i=0) const
G4GLOB_DLL std::ostream G4cout
G4double GetUserElapsed() const
Definition: G4Timer.cc:130
G4PrimaryVertex * GetPrimaryVertex(G4int i=0) const
Definition: G4Event.hh:167
G4double GetRealElapsed() const
Definition: G4Timer.cc:107
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void Stop()
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
double y() const
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61
Hep3Vector cross(const Hep3Vector &) const
G4HCofThisEvent * GetHCofThisEvent() const
Definition: G4Event.hh:185
double G4double
Definition: G4Types.hh:76
double mag() const

Here is the call graph for this function:


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