File indexing completed on 2025-01-18 09:17:11
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
0030
0031
0032
0033
0034
0035 #include "XraySPOSteppingAction.hh"
0036 #include "XraySPODetectorConstruction.hh"
0037 #include "XraySPOHistoManager.hh"
0038 #include "G4SteppingManager.hh"
0039 #include "G4RunManager.hh"
0040 #include "Randomize.hh"
0041
0042
0043
0044 void XraySPOSteppingAction::UserSteppingAction(const G4Step* step)
0045 {
0046 G4int parentID = step->GetTrack()->GetParentID();
0047
0048
0049 if (parentID == 0)
0050 {
0051 G4RunManager *rm = G4RunManager::GetRunManager();
0052 G4int eventID = rm->GetCurrentEvent()->GetEventID();
0053
0054 G4double len_fraction = 0.1;
0055 G4bool proceed = false;
0056
0057 if (fPrevEventID != eventID)
0058 {
0059 fNumReflections = -2;
0060 }
0061
0062 G4ThreeVector direction = step->GetTrack()->GetMomentumDirection();
0063 auto theta = (G4double)direction.theta();
0064 auto phi = (G4double)direction.phi();
0065
0066 G4double x = step->GetPreStepPoint()->GetPosition().x();
0067 G4double y = step->GetPreStepPoint()->GetPosition().y();
0068 G4double z = step->GetPreStepPoint()->GetPosition().z();
0069
0070 const G4String& vol_name = step->GetTrack()->GetVolume()->GetName();
0071 G4int trackid = step->GetTrack()->GetTrackID();
0072
0073 G4StepPoint * postStep = step->GetPostStepPoint();
0074 G4StepPoint* post_point = step->GetPostStepPoint();
0075
0076 const G4String& proc_name = postStep->GetProcessDefinedStep()->GetProcessName();
0077
0078
0079 G4String s = "";
0080 G4String& post_phys_name = s;
0081 if (fPrevEventID == eventID && step->GetStepLength() >= len_fraction)
0082 {
0083 if (post_point != nullptr)
0084 {
0085 G4VPhysicalVolume* post_phys = post_point->GetPhysicalVolume();
0086 if (post_phys != nullptr)
0087 {
0088 post_phys_name = post_phys->GetName();
0089 }
0090
0091 if (post_phys_name == "cone_1" || post_phys_name == "cone_2" || post_phys_name == "cone_3" || post_phys_name == "cone_4" || post_phys_name == "cone_5" || post_phys_name == "cone_6" || post_phys_name == "cone_7" || post_phys_name == "cone_8")
0092 {
0093 if (fNumReflections < 0)
0094 {
0095 fNumReflections = 1;
0096 }
0097 else
0098 {
0099 fNumReflections += 1;
0100 }
0101 proceed = true;
0102 }
0103 }
0104 }
0105
0106 if (vol_name == "pDummyEntrance" || vol_name == "pDummyExit" || vol_name == "pDummySphere")
0107 {
0108 proceed = true;
0109 }
0110
0111 if (proceed)
0112 {
0113
0114 G4AnalysisManager* analysisManager = G4AnalysisManager::Instance();
0115 analysisManager->FillNtupleIColumn(0, eventID);
0116 analysisManager->FillNtupleSColumn(1, vol_name);
0117 analysisManager->FillNtupleIColumn(2, trackid);
0118 analysisManager->FillNtupleDColumn(3, x);
0119 analysisManager->FillNtupleDColumn(4, y);
0120 analysisManager->FillNtupleDColumn(5, z);
0121 analysisManager->FillNtupleDColumn(6, theta);
0122 analysisManager->FillNtupleDColumn(7, phi);
0123 analysisManager->FillNtupleSColumn(8, proc_name);
0124 analysisManager->FillNtupleIColumn(9, parentID);
0125 analysisManager->FillNtupleIColumn(10, fNumReflections);
0126 analysisManager->AddNtupleRow();
0127 }
0128 fPrevx = x;
0129 fPrevy = y;
0130 fPrevz = z;
0131 fPrevEventID = eventID;
0132 }
0133 }
0134
0135