Warning, file /geant4/examples/extended/exoticphysics/saxs/src/SAXSSteppingAction.cc was not indexed
or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029 #ifdef G4MULTITHREADED
0030 # include "G4MTRunManager.hh"
0031 #else
0032 # include "G4RunManager.hh"
0033 #endif
0034
0035 #include "SAXSDetectorConstruction.hh"
0036 #include "SAXSEventAction.hh"
0037 #include "SAXSSteppingAction.hh"
0038
0039 #include "G4AnalysisManager.hh"
0040 #include "G4Event.hh"
0041 #include "G4EventManager.hh"
0042 #include "G4HCofThisEvent.hh"
0043 #include "G4OpticalPhoton.hh"
0044 #include "G4SDManager.hh"
0045 #include "G4Step.hh"
0046 #include "G4SystemOfUnits.hh"
0047 #include "G4Track.hh"
0048 #include "G4Trajectory.hh"
0049 #include "G4TrajectoryContainer.hh"
0050 #include "G4UImanager.hh"
0051 #include "G4VHitsCollection.hh"
0052 #include "G4VVisManager.hh"
0053 #include "G4ios.hh"
0054
0055 #include <cmath>
0056
0057
0058
0059 SAXSSteppingAction::SAXSSteppingAction(SAXSEventAction* eventAction)
0060 : G4UserSteppingAction(), fEventAction(eventAction), fPhantom(0), fEventNumber(-1), fNSe(0)
0061 {}
0062
0063
0064
0065 SAXSSteppingAction::~SAXSSteppingAction() {}
0066
0067
0068
0069 void SAXSSteppingAction::UserSteppingAction(const G4Step* step)
0070 {
0071
0072 const SAXSDetectorConstruction* detectorConstruction =
0073 static_cast<const SAXSDetectorConstruction*>(
0074 G4RunManager::GetRunManager()->GetUserDetectorConstruction());
0075
0076
0077 fPhantom = detectorConstruction->GetPhantom();
0078
0079
0080 G4StepPoint* preStepPoint = step->GetPreStepPoint();
0081 G4StepPoint* postStepPoint = step->GetPostStepPoint();
0082
0083
0084 G4LogicalVolume* volume = preStepPoint->GetTouchableHandle()->GetVolume()->GetLogicalVolume();
0085
0086
0087 G4Track* track = step->GetTrack();
0088 G4int ID = track->GetTrackID();
0089 G4String partName = track->GetDefinition()->GetParticleName();
0090
0091
0092 G4double PrimaryWeight = 1.;
0093 if (partName == "gamma" && ID == 1) {
0094 PrimaryWeight = track->GetWeight();
0095 fEventAction->UpdateEventWeight(PrimaryWeight);
0096 }
0097
0098
0099 if (volume != fPhantom || partName != "gamma" || ID != 1) return;
0100
0101
0102 G4int eventNumber = G4RunManager::GetRunManager()->GetCurrentEvent()->GetEventID();
0103 if (eventNumber != fEventNumber) {
0104 fEventNumber = eventNumber;
0105 fNSe = 0;
0106 }
0107
0108
0109 G4String Process;
0110 if (postStepPoint->GetProcessDefinedStep() != NULL) {
0111 Process = postStepPoint->GetProcessDefinedStep()->GetProcessName();
0112 }
0113 else {
0114 Process = "UserLimit";
0115 }
0116
0117 G4int ProcIndex = -1;
0118 if (Process == "Transportation" || Process == "CoupledTransportation" || Process == "UserLimit") {
0119 ProcIndex = 0;
0120 }
0121 if (Process == "Rayl" || Process == "biasWrapper(Rayl)") {
0122 ProcIndex = 1;
0123 fNSe++;
0124 fEventAction->AddNRi();
0125 }
0126 if (Process == "compt" || Process == "biasWrapper(compt)") {
0127 ProcIndex = 2;
0128 fNSe++;
0129 fEventAction->AddNCi();
0130 }
0131 if (Process == "phot" || Process == "biasWrapper(phot)") ProcIndex = 3;
0132 if (Process == "conv" || Process == "biasWrapper(conv)") ProcIndex = 4;
0133 if (Process == "photonNuclear" || Process == "biasWrapper(photonNuclear)") {
0134 ProcIndex = 5;
0135 }
0136 if (Process == "PowderDiffraction" || Process == "biasWrapper(PowderDiffraction)") {
0137 ProcIndex = 6;
0138 fNSe++;
0139 fEventAction->AddNDi();
0140 }
0141
0142
0143
0144 G4ThreeVector mom1 = preStepPoint->GetMomentumDirection();
0145 G4ThreeVector mom2 = postStepPoint->GetMomentumDirection();
0146 G4double theta = mom1.angle(mom2);
0147
0148
0149 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
0150 if (theta > 1e-6) {
0151 analysisManager->FillNtupleIColumn(1, 0, ProcIndex);
0152 analysisManager->FillNtupleDColumn(1, 1, preStepPoint->GetKineticEnergy() / CLHEP::keV);
0153 analysisManager->FillNtupleDColumn(1, 2, theta / CLHEP::deg);
0154 analysisManager->FillNtupleDColumn(1, 3, track->GetWeight());
0155 analysisManager->AddNtupleRow(1);
0156 }
0157 }
0158
0159