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

#include <ExExChEventAction.hh>

Inheritance diagram for ExExChEventAction:
Collaboration diagram for ExExChEventAction:

Public Member Functions

 ExExChEventAction ()
 
virtual ~ExExChEventAction ()
 
virtual void BeginOfEventAction (const G4Event *)
 
virtual void EndOfEventAction (const G4Event *)
 
void SetVerbose (G4int val)
 
G4int GetVerbose () const
 
- 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 41 of file ExExChEventAction.hh.

Constructor & Destructor Documentation

ExExChEventAction::ExExChEventAction ( )

Definition at line 51 of file ExExChEventAction.cc.

52 {
53  fSD_ID = -1;
54  fVerboseLevel = 0;
55 }
ExExChEventAction::~ExExChEventAction ( )
virtual

Definition at line 59 of file ExExChEventAction.cc.

59  {
60 }

Member Function Documentation

void ExExChEventAction::BeginOfEventAction ( const G4Event )
virtual

Reimplemented from G4UserEventAction.

Definition at line 64 of file ExExChEventAction.cc.

64  {
65 }
void ExExChEventAction::EndOfEventAction ( const G4Event evt)
virtual

Reimplemented from G4UserEventAction.

Definition at line 69 of file ExExChEventAction.cc.

70 {
72 
73  if(fSD_ID == -1) {
74  G4String sdName;
75  if(SDman->FindSensitiveDetector(sdName="telescope",0)){
76  fSD_ID = SDman->GetCollectionID(sdName="telescope/collection");
77  }
78  }
79 
81 
82  G4HCofThisEvent * HCE = evt->GetHCofThisEvent();
83 
84  if(HCE)
85  {
86  if(fSD_ID != -1){
87  G4VHitsCollection* aHCSD = HCE->GetHC(fSD_ID);
89  }
90  }
91 
92  G4ThreeVector SSDposition[3];
93  SSDposition[0]= G4ThreeVector(0.,0.,0.);
94  SSDposition[1]= G4ThreeVector(0.,0.,0.);
95  SSDposition[2]= G4ThreeVector(0.,0.,0.);
96 
97  if(fSD)
98  {
99  int bTotalHits = 0;
100 
101  int n_hit_sd = fSD->entries();
102  for(int i2=0;i2<3;i2++){
103  for(int i1=0;i1<n_hit_sd;i1++)
104  {
105  ExExChSensitiveDetectorHit* aHit = (*fSD)[i1];
106  if(aHit->GetLayerID()==i2) {
107  SSDposition[i2] = aHit->GetWorldPos();
108  bTotalHits++;
109  }
110  }
111  }
112 
113  if(bTotalHits == 3){
114  double fAngXin = (SSDposition[1].x() - SSDposition[0].x());
115  fAngXin /= (SSDposition[1].z() - SSDposition[0].z());
116  double fAngYin = (SSDposition[1].y() - SSDposition[0].y());
117  fAngYin /= (SSDposition[1].z() - SSDposition[0].z());
118  double fPosXin = SSDposition[1].x();
119  double fPosYin = SSDposition[1].y();
120  double fAngXout = (SSDposition[2].x() - SSDposition[1].x());
121  fAngXout /= (SSDposition[2].z() - SSDposition[1].z());
122  double fAngYout = (SSDposition[2].y() - SSDposition[1].y());
123  fAngYout /= (SSDposition[2].z() - SSDposition[1].z());
124 
125  G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
126 
127  analysisManager->FillNtupleDColumn(0, fAngXin * 1.E6 * CLHEP::rad);
128  analysisManager->FillNtupleDColumn(1, fAngYin * 1.E6 * CLHEP::rad);
129  analysisManager->FillNtupleDColumn(2, fPosXin / CLHEP::mm);
130  analysisManager->FillNtupleDColumn(3, fPosYin / CLHEP::mm);
131  analysisManager->FillNtupleDColumn(4, fAngXout * 1.E6 * CLHEP::rad);
132  analysisManager->FillNtupleDColumn(5, fAngYout * 1.E6 * CLHEP::rad);
133  analysisManager->AddNtupleRow();
134  }
135  }
136 }
G4VHitsCollection * GetHC(G4int i)
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:135
CLHEP::Hep3Vector G4ThreeVector
double x() const
static constexpr double rad
G4int entries() const
double z() const
G4bool FillNtupleDColumn(G4int id, G4double value)
static constexpr double mm
Definition: SystemOfUnits.h:95
G4VSensitiveDetector * FindSensitiveDetector(G4String dName, G4bool warning=true)
Definition: G4SDManager.cc:128
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
double y() const
G4HCofThisEvent * GetHCofThisEvent() const
Definition: G4Event.hh:185

Here is the call graph for this function:

G4int ExExChEventAction::GetVerbose ( ) const
inline

Definition at line 56 of file ExExChEventAction.hh.

56 { return fVerboseLevel; }
void ExExChEventAction::SetVerbose ( G4int  val)
inline

Definition at line 55 of file ExExChEventAction.hh.

55 { fVerboseLevel = val; }

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