Loading...
Searching...
No Matches
SAXSEventAction.cc
Go to the documentation of this file.
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/// \file SAXSEventAction.cc
27/// \brief Implementation of the SAXSEventAction class
28//
29//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
30
31#ifdef G4MULTITHREADED
32#include "G4MTRunManager.hh"
33#else
34#include "G4RunManager.hh"
35#endif
36
37#include "G4Event.hh"
38#include "G4EventManager.hh"
39#include "G4HCofThisEvent.hh"
40#include "G4VHitsCollection.hh"
41#include "G4TrajectoryContainer.hh"
42#include "G4Trajectory.hh"
43#include "G4VVisManager.hh"
44#include "G4SDManager.hh"
45#include "G4UImanager.hh"
46#include "G4ios.hh"
47#include "G4SystemOfUnits.hh"
48#include "G4AnalysisManager.hh"
49
50#include "SAXSEventAction.hh"
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54
57 fSensitiveDetector_ID(-1),
58 fVerboseLevel(0),
59 fNRi(0),
60 fNCi(0),
61 fNDi(0),
62 fEventWeight(0.)
63{}
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
70
72{
73 fNRi = 0;
74 fNCi = 0;
75 fNDi = 0;
76 fEventWeight = 1.;
77}
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80
82{
83 //instantiating The Sensitive Detector Manager
84 G4SDManager* SDman = G4SDManager::GetSDMpointer();
85
86 //Hit Detection System
87 if (fSensitiveDetector_ID==-1) {
88 G4String SensitiveDetectorName;
89 if (SDman->FindSensitiveDetector(SensitiveDetectorName="det",0)) {
91 SDman->GetCollectionID(SensitiveDetectorName="det/collection");
92 }
93 }
94
95 SensitiveDetectorHitsCollection* fSensitiveDetectorHC = 0;
96 G4HCofThisEvent* HCE = aEvent->GetHCofThisEvent();
97
98 if (HCE) {
99 if (fSensitiveDetector_ID != -1) {
100 G4VHitsCollection* aHC = HCE->GetHC(fSensitiveDetector_ID);
101 fSensitiveDetectorHC = (SensitiveDetectorHitsCollection*)(aHC);
102 }
103 }
104
105 //instantiating The Analysis Manager
106 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
107
108 //get the Event Number
109 G4int eventNumber =
110 G4RunManager::GetRunManager()->GetCurrentEvent()->GetEventID();
111
112 //filling the SD scoring ntuple
113 if (fSensitiveDetectorHC) {
114 size_t vNumberOfHit = fSensitiveDetectorHC->entries();
115 for (size_t i=0; i<vNumberOfHit; i++) {
116 SAXSSensitiveDetectorHit* aHit = (*fSensitiveDetectorHC)[i];
117 analysisManager->FillNtupleDColumn(0,0,aHit->GetEnergy()/CLHEP::keV);
118 analysisManager->FillNtupleDColumn(0,1,aHit->GetPos().x()/CLHEP::mm);
119 analysisManager->FillNtupleDColumn(0,2,aHit->GetPos().y()/CLHEP::mm);
120 analysisManager->FillNtupleDColumn(0,3,aHit->GetPos().z()/CLHEP::mm);
121 analysisManager->FillNtupleDColumn(0,4,aHit->GetMom().x());
122 analysisManager->FillNtupleDColumn(0,5,aHit->GetMom().y());
123 analysisManager->FillNtupleDColumn(0,6,aHit->GetMom().z());
124 analysisManager->FillNtupleDColumn(0,7,aHit->GetTime()/CLHEP::ns);
125 analysisManager->FillNtupleIColumn(0,8,aHit->GetType());
126 analysisManager->FillNtupleIColumn(0,9,aHit->GetTrackID());
127 analysisManager->FillNtupleIColumn(0,10,fNRi);
128 analysisManager->FillNtupleIColumn(0,11,fNCi);
129 analysisManager->FillNtupleIColumn(0,12,fNDi);
130 analysisManager->FillNtupleIColumn(0,13,eventNumber);
131 analysisManager->FillNtupleDColumn(0,14,aHit->GetWeight());
132 analysisManager->AddNtupleRow(0);
133 }
134 }
135
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
139
Implementation of the SAXSEventAction class.
Definition of the SAXSSensitiveDetectorHit class.
G4THitsCollection< SAXSSensitiveDetectorHit > SensitiveDetectorHitsCollection
virtual ~SAXSEventAction()
virtual void EndOfEventAction(const G4Event *)
virtual void BeginOfEventAction(const G4Event *)

Applications | User Support | Publications | Collaboration