Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
SAXSSteppingAction Class Reference

Stepping Action. More...

#include <Doxymodules_exoticphysics.h>

Inheritance diagram for SAXSSteppingAction:
G4UserSteppingAction

Public Member Functions

 SAXSSteppingAction (SAXSEventAction *)
 
virtual ~SAXSSteppingAction ()
 
virtual void UserSteppingAction (const G4Step *)
 

Private Attributes

SAXSEventActionfEventAction
 
G4LogicalVolumefPhantom
 
G4int fEventNumber
 
G4int fNSe
 

Detailed Description

Stepping Action.

It is used to retrive information about every interaction of the photons occuring iside the phantom, and in particular, scattering events. The total number of scattering events (subdivided by tyoe) is passed to the EventAction class and stored at the end of the event.

Definition at line 75 of file Doxymodules_exoticphysics.h.

Constructor & Destructor Documentation

◆ SAXSSteppingAction()

SAXSSteppingAction::SAXSSteppingAction ( SAXSEventAction eventAction)

Definition at line 62 of file SAXSSteppingAction.cc.

62 :
64 fEventAction(eventAction),
65 fPhantom(0),
66 fEventNumber(-1),
67 fNSe(0)
68{}
G4LogicalVolume * fPhantom
SAXSEventAction * fEventAction

◆ ~SAXSSteppingAction()

SAXSSteppingAction::~SAXSSteppingAction ( )
virtual

Definition at line 72 of file SAXSSteppingAction.cc.

72{}

Member Function Documentation

◆ UserSteppingAction()

void SAXSSteppingAction::UserSteppingAction ( const G4Step step)
virtual

Definition at line 76 of file SAXSSteppingAction.cc.

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}
G4LogicalVolume * GetPhantom() const
void UpdateEventWeight(G4double ew)

Member Data Documentation

◆ fEventAction

SAXSEventAction* SAXSSteppingAction::fEventAction
private

Definition at line 58 of file SAXSSteppingAction.hh.

◆ fPhantom

G4LogicalVolume* SAXSSteppingAction::fPhantom
private

Definition at line 60 of file SAXSSteppingAction.hh.

◆ fEventNumber

G4int SAXSSteppingAction::fEventNumber
private

Definition at line 62 of file SAXSSteppingAction.hh.

◆ fNSe

G4int SAXSSteppingAction::fNSe
private

Definition at line 63 of file SAXSSteppingAction.hh.


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

Applications | User Support | Publications | Collaboration