File indexing completed on 2025-02-23 09:20:36
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
0036 #include "XrayTESdetSteppingAction.hh"
0037 #include "G4AnalysisManager.hh"
0038
0039 #include "G4SteppingManager.hh"
0040 #include "G4RunManager.hh"
0041 #include "G4SystemOfUnits.hh"
0042 #include "G4PhysicalConstants.hh"
0043
0044
0045
0046 G4String clean_name (const G4String& vol_name)
0047 {
0048 G4String cleaned_name = "";
0049
0050
0051 if (((std::string)vol_name).find("Bipxl") != std::string::npos)
0052 cleaned_name = "Bipxl";
0053 if (((std::string)vol_name).find("membrane") != std::string::npos)
0054 cleaned_name = "membrane";
0055 if (((std::string)vol_name).find("gridpiece") != std::string::npos)
0056 cleaned_name = "gridpiece";
0057 if (((std::string)vol_name).find("trapezoid") != std::string::npos)
0058 cleaned_name = "Mesh";
0059 if (cleaned_name == "") cleaned_name = vol_name;
0060 return cleaned_name;
0061 }
0062
0063
0064
0065 void XrayTESdetSteppingAction::UserSteppingAction(const G4Step* step)
0066 {
0067
0068 G4RunManager *rm = G4RunManager::GetRunManager();
0069 auto analysisManager = G4AnalysisManager::Instance();
0070
0071 G4int eventID = rm->GetCurrentEvent()->GetEventID();
0072 G4Track* track = step->GetTrack();
0073 G4StepPoint* pre_step_point = step->GetPreStepPoint();
0074 G4int parentID = track->GetParentID();
0075 G4int trackID = track->GetTrackID();
0076 const G4String& vol_name = track->GetVolume()->GetName();
0077 G4String mother_name = "";
0078 G4double init_kinetic_energy = 0;
0079 G4double kinetic_energy = track->GetKineticEnergy();
0080 G4double pre_kinetic_energy = pre_step_point->GetKineticEnergy();
0081 const G4String& particle_name = track->GetParticleDefinition()->GetParticleName();
0082 const G4VProcess* pre_step_proc = pre_step_point->GetProcessDefinedStep();
0083 const G4VProcess* post_step_proc = pre_step_point->GetProcessDefinedStep();
0084
0085 bool proceed = false;
0086 G4String pre_step_name_s = "Undefined";
0087 G4String post_step_name_s = "Undefined";
0088 G4String creator_process_name = "Undefined";
0089 G4String init_creator_process = "Undefined";
0090 const G4String& pre_step_name = pre_step_name_s;
0091 const G4String& post_step_name = post_step_name_s;
0092
0093 if (track->GetCreatorProcess() != nullptr)
0094 {
0095 creator_process_name = track->GetCreatorProcess()->GetProcessName();
0096 }
0097
0098 if (pre_step_proc != nullptr)
0099 {
0100 pre_step_name_s = pre_step_proc->GetProcessName();
0101 }
0102
0103 if (post_step_proc != nullptr)
0104 {
0105 post_step_name_s = post_step_proc->GetProcessName();
0106 }
0107
0108
0109 if (eventID != fPrev_eventID)
0110 {
0111 fInit_energy.clear();
0112 fCreator_proc.clear();
0113 }
0114
0115 G4int current_step_number = -1;
0116 current_step_number = track->GetCurrentStepNumber();
0117 if (current_step_number != -1)
0118 {
0119 if (current_step_number == 1 && parentID == 0)
0120 {
0121 proceed = true;
0122 fInit_energy.insert(std::pair<G4int, G4double>(trackID, pre_kinetic_energy));
0123 fCreator_proc.insert(std::pair<G4int, G4String>(trackID, creator_process_name));
0124 pre_step_name_s = "InitStep";
0125 }
0126 }
0127
0128
0129 if (trackID > 1)
0130 {
0131
0132 if (fInit_energy.count((int)trackID) == 0)
0133 {
0134 fInit_energy.insert(std::pair<G4int, G4double>(trackID, pre_kinetic_energy));
0135 fCreator_proc.insert(std::pair<G4int, G4String>(trackID, creator_process_name));
0136 }
0137 }
0138
0139
0140 if (((std::string)vol_name).find("expHall") == std::string::npos)
0141 {
0142 if (track->GetVolume() != nullptr)
0143 {
0144 if (track->GetVolume()->GetMotherLogical() != nullptr)
0145 {
0146 mother_name = track->GetVolume()->GetMotherLogical()->GetName();
0147 }
0148 }
0149 if (((std::string)mother_name).find("Detector_log") != std::string::npos || ((std::string)mother_name).find("element_log") != std::string::npos)
0150 {
0151 proceed = true;
0152 }
0153 }
0154
0155 init_kinetic_energy = fInit_energy[(int)trackID];
0156 init_creator_process = fCreator_proc[trackID];
0157
0158 G4String cleaned_name = "";
0159 cleaned_name = clean_name(vol_name);
0160
0161
0162 if (proceed)
0163 {
0164 G4double x = pre_step_point->GetPosition().x();
0165 G4double y = pre_step_point->GetPosition().y();
0166 G4double z = pre_step_point->GetPosition().z();
0167 G4ThreeVector direction = pre_step_point->GetMomentum();
0168 G4double theta = direction.theta();
0169 G4double phi = direction.phi();
0170 G4int pixel_number = track->GetVolume()->GetCopyNo();
0171 G4int def_replica_number = -999;
0172 G4int step_number = track->GetCurrentStepNumber();
0173 G4double step_energy_dep = step->GetTotalEnergyDeposit();
0174
0175
0176 analysisManager->FillNtupleIColumn(0, eventID);
0177 analysisManager->FillNtupleSColumn(1, cleaned_name);
0178 analysisManager->FillNtupleIColumn(2, trackID);
0179 analysisManager->FillNtupleDColumn(3, x);
0180 analysisManager->FillNtupleDColumn(4, y);
0181 analysisManager->FillNtupleDColumn(5, z);
0182 analysisManager->FillNtupleDColumn(6, theta);
0183 analysisManager->FillNtupleDColumn(7, phi);
0184 analysisManager->FillNtupleIColumn(8, parentID);
0185 if (((std::string)vol_name).find("Bipxl") != std::string::npos)
0186 {
0187 analysisManager->FillNtupleIColumn(9, pixel_number);
0188 }
0189 else
0190 {
0191 analysisManager->FillNtupleIColumn(9, def_replica_number);
0192 }
0193 analysisManager->FillNtupleDColumn(10, step_energy_dep);
0194 analysisManager->FillNtupleIColumn(11, step_number);
0195 analysisManager->FillNtupleDColumn(12, init_kinetic_energy);
0196 analysisManager->FillNtupleDColumn(13, kinetic_energy);
0197 analysisManager->FillNtupleSColumn(14, particle_name);
0198 analysisManager->FillNtupleSColumn(15, pre_step_name);
0199 analysisManager->FillNtupleSColumn(16, post_step_name);
0200 analysisManager->FillNtupleSColumn(17, init_creator_process);
0201 analysisManager->AddNtupleRow(0);
0202 }
0203 fPrev_eventID = eventID;
0204 fPrev_trackID = trackID;
0205 }
0206
0207