Loading...
Searching...
No Matches
SAXSSteppingAction.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 SAXSSteppingAction.cc
27/// \brief Definition of the SAXSSteppingAction 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 "SAXSSteppingAction.hh"
51#include "SAXSEventAction.hh"
53
54#include "G4Step.hh"
55#include "G4Track.hh"
56#include "G4OpticalPhoton.hh"
57
58#include <cmath>
59
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
61
64 fEventAction(eventAction),
65 fPhantom(0),
66 fEventNumber(-1),
67 fNSe(0)
68{}
69
70//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
71
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
75
77{
78 //DetectorConstruction instance
79 const SAXSDetectorConstruction* detectorConstruction =
80 static_cast<const SAXSDetectorConstruction*>
81 (G4RunManager::GetRunManager()->GetUserDetectorConstruction());
82
83 //get Phantom volume
84 fPhantom = detectorConstruction->GetPhantom();
85
86 //get pre and post step points
87 G4StepPoint* preStepPoint = step->GetPreStepPoint();
88 G4StepPoint* postStepPoint = step->GetPostStepPoint();
89
90 //get the volume of the current step
91 G4LogicalVolume* volume = preStepPoint->GetTouchableHandle()
92 ->GetVolume()->GetLogicalVolume();
93
94 //get track and particle name
95 G4Track* track = step->GetTrack();
96 G4int ID = track->GetTrackID();
97 G4String partName = track->GetDefinition()->GetParticleName();
98
99 //update the primary particle (event) weight
100 G4double PrimaryWeight = 1.;
101 if (partName == "gamma" && ID==1) {
102 PrimaryWeight = track->GetWeight();
103 fEventAction->UpdateEventWeight(PrimaryWeight);
104 }
105
106 //if the particle is not a primary photon stepping in the phantom, exit!
107 if (volume != fPhantom || partName != "gamma" || ID!=1) return;
108
109 //if it is a new event, reset variable values to 0
110 G4int eventNumber = G4RunManager::GetRunManager()
111 ->GetCurrentEvent()->GetEventID();
112 if (eventNumber != fEventNumber) {
113 fEventNumber = eventNumber;
114 fNSe = 0;
115 }
116
117 //identify the process that occurred
118 G4String Process;
119 if (postStepPoint->GetProcessDefinedStep() != NULL) {
120 Process = postStepPoint->GetProcessDefinedStep()->GetProcessName();
121 } else {
122 Process = "UserLimit";
123 }
124 //G4cout << "Process: " << Process << G4endl;
125 G4int ProcIndex = -1;
126 if (Process == "Transportation" || Process == "CoupledTransportation"
127 || Process == "UserLimit") {
128 ProcIndex = 0;
129 }
130 if (Process == "Rayl" || Process == "biasWrapper(Rayl)") {
131 ProcIndex = 1;
132 fNSe++;
134 }
135 if (Process == "compt" || Process == "biasWrapper(compt)") {
136 ProcIndex = 2;
137 fNSe++;
139 }
140 if (Process == "phot" || Process == "biasWrapper(phot)") ProcIndex = 3;
141 if (Process == "conv" || Process == "biasWrapper(conv)") ProcIndex = 4;
142 if (Process == "photonNuclear" ||
143 Process == "biasWrapper(photonNuclear)") {
144 ProcIndex = 5;
145 }
146 if (Process == "PowderDiffraction" ||
147 Process == "biasWrapper(PowderDiffraction)") {
148 ProcIndex = 6;
149 fNSe++;
151 }
152 //G4cout << "ProcIndex: " << ProcIndex << G4endl;
153
154 //calculate the scattering angle
155 G4ThreeVector mom1 = preStepPoint->GetMomentumDirection();
156 G4ThreeVector mom2 = postStepPoint->GetMomentumDirection();
157 G4double theta = mom1.angle(mom2);
158
159 //fill the Scattering ntuple
160 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
161 if (theta > 1e-6) { //record only if the angle is not 0
162 analysisManager->FillNtupleIColumn(1,0,ProcIndex);
163 analysisManager->FillNtupleDColumn(1,1,
164 preStepPoint->GetKineticEnergy()/CLHEP::keV);
165 analysisManager->FillNtupleDColumn(1,2,theta/CLHEP::deg);
166 analysisManager->FillNtupleDColumn(1,3,track->GetWeight());
167 analysisManager->AddNtupleRow(1);
168 }
169}
170
171//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
172
Implementation of the SAXSDetectorConstruction class.
Implementation of the SAXSEventAction class.
Definition of the SAXSSteppingAction class.
G4LogicalVolume * GetPhantom() const
void UpdateEventWeight(G4double ew)
virtual void UserSteppingAction(const G4Step *)
G4LogicalVolume * fPhantom
SAXSSteppingAction(SAXSEventAction *)
SAXSEventAction * fEventAction

Applications | User Support | Publications | Collaboration