9#include <simulation/kernel/SteppingAction.h>
10#include <simulation/kernel/EventAction.h>
11#include <simulation/kernel/RunManager.h>
12#include <simulation/kernel/UserInfo.h>
13#include <framework/core/Environment.h>
14#include <framework/logging/Logger.h>
15#include <framework/gearbox/Unit.h>
18#include <G4UnitsTable.hh>
19#include <G4VProcess.hh>
20#include <G4VProcess.hh>
25using namespace Simulation;
27SteppingAction::SteppingAction()
30 m_maxNumberSteps = 100000;
33 UserSteppingAction(aStep);
39SteppingAction::~SteppingAction()
45void SteppingAction::UserSteppingAction(
const G4Step* step)
47 G4Track* track = step->GetTrack();
52 if (track->GetVolume() == NULL) {
53 B2WARNING(
"SteppingAction: Track in NULL volume, terminating!\n"
54 <<
"step_no=" << track->GetCurrentStepNumber() <<
" type=" << track->GetDefinition()->GetParticleName()
55 <<
"\n position=" << G4BestUnit(track->GetPosition(),
"Length") <<
" momentum=" << G4BestUnit(track->GetMomentum(),
"Energy"));
56 track->SetTrackStatus(fStopAndKill);
63 for (
auto& rAbsorber : m_absorbers) {
64 const G4ThreeVector stepPrePos = step->GetPreStepPoint()->GetPosition() / CLHEP::mm *
Unit::mm;
65 const G4ThreeVector stepPostPos = step->GetPostStepPoint()->GetPosition() / CLHEP::mm *
Unit::mm;
66 if (stepPrePos.perp() < (rAbsorber *
Unit::cm) && stepPostPos.perp() > (rAbsorber *
Unit::cm)) {
70 track->SetTrackStatus(fStopAndKill);
76 if (m_writeSimSteps) {
77 writeVREventStep(step, track);
83 if (track->GetCurrentStepNumber() > m_maxNumberSteps) {
84 B2WARNING(
"SteppingAction: Too many steps for this track, terminating!\n"
85 <<
"step_no=" << track->GetCurrentStepNumber() <<
"type=" << track->GetDefinition()->GetParticleName()
86 <<
"\n position=" << G4BestUnit(track->GetPosition(),
"Length") <<
" momentum=" << G4BestUnit(track->GetMomentum(),
"Energy"));
87 track->SetTrackStatus(fStopAndKill);
94 if (m_storeTrajectories) {
96 if (info && info->getTrajectory()) {
98 if (trajectory.
empty()) {
99 const G4ThreeVector stepPos = step->GetPreStepPoint()->GetPosition() / CLHEP::mm *
Unit::mm;
100 const G4ThreeVector stepMom = step->GetPreStepPoint()->GetMomentum() / CLHEP::MeV *
Unit::MeV;
102 stepPos.x(), stepPos.y(), stepPos.z(),
103 stepMom.x(), stepMom.y(), stepMom.z()
106 const G4ThreeVector stepPos = step->GetPostStepPoint()->GetPosition() / CLHEP::mm *
Unit::mm;
107 const G4ThreeVector stepMom = step->GetPostStepPoint()->GetMomentum() / CLHEP::MeV *
Unit::MeV;
109 stepPos.x(), stepPos.y(), stepPos.z(),
110 stepMom.x(), stepMom.y(), stepMom.z()
117void SteppingAction::writeVREventStep(
const G4Step* step,
const G4Track* track)
120 RunManager& runManager = RunManager::Instance();
122 if (eventAction ==
nullptr)
125 if (output ==
nullptr)
127 if (!output->is_open())
130 if (abs(track->GetDefinition()->GetPDGEncoding()) > 1000020040)
133 G4StepPoint* postStepPoint = step->GetPostStepPoint();
134 if (postStepPoint->GetGlobalTime() > 100.0)
137 G4StepPoint* preStepPoint = step->GetPreStepPoint();
138 double KE = preStepPoint->GetTotalEnergy() - track->GetDefinition()->GetPDGMass();
139 if ((KE < 0.0005) && (track->GetDefinition()->GetParticleName() !=
"opticalphoton"))
143 G4String pVolName = track->GetVolume()->GetName();
144 G4String sensitiveDetectorName =
"";
145 if (pVolName.compare(0, 4,
"PXD.") == 0) {
146 if (pVolName.find(
".Active") != std::string::npos) { sensitiveDetectorName =
"PXD"; }
147 }
else if (pVolName.compare(0, 4,
"SVD.") == 0) {
148 if (pVolName.find(
".Active") != std::string::npos) { sensitiveDetectorName =
"SVD"; }
149 }
else if (pVolName.compare(0, 20,
"physicalSD_CDCLayer_") == 0) {
150 sensitiveDetectorName =
"CDC";
151 }
else if (pVolName.compare(0, 19,
"TOP.moduleSensitive") == 0) {
152 sensitiveDetectorName =
"TOP";
153 }
else if (pVolName.compare(0, 23,
"av_1_impr_1_cuttest_pv_") == 0) {
154 sensitiveDetectorName =
"TOP";
155 }
else if (pVolName.compare(0, 12,
"moduleWindow") == 0) {
156 sensitiveDetectorName =
"ARICH";
157 }
else if (pVolName.compare(0, 25,
"ARICH.AerogelSupportPlate") == 0) {
158 sensitiveDetectorName =
"ARICH";
159 }
else if (pVolName.compare(0, 25,
"eclBarrelCrystalPhysical_") == 0) {
160 sensitiveDetectorName =
"ECL";
161 }
else if (pVolName.compare(0, 22,
"eclFwdCrystalPhysical_") == 0) {
162 sensitiveDetectorName =
"ECL";
163 }
else if (pVolName.compare(0, 22,
"eclBwdCrystalPhysical_") == 0) {
164 sensitiveDetectorName =
"ECL";
165 }
else if (pVolName.compare(0, 20,
"BKLM.ScintActiveType") == 0) {
166 sensitiveDetectorName =
"BKLM";
167 }
else if (pVolName.compare(0, 10,
"BKLM.Layer") == 0) {
168 if (pVolName.find(
"GasPhysical") != std::string::npos) {
169 sensitiveDetectorName =
"BKLM";
171 }
else if (pVolName.compare(0, 15,
"StripSensitive_") == 0) {
172 sensitiveDetectorName =
"EKLM";
180 (*output) << std::fixed << std::setprecision(4)
181 << track->GetTrackID() <<
","
182 << track->GetParentID() <<
","
183 << track->GetDefinition()->GetParticleName() <<
","
184 << track->GetDefinition()->GetPDGMass() <<
","
185 << int(track->GetDefinition()->GetPDGCharge()) <<
","
186 << track->GetCurrentStepNumber() <<
","
187 << track->GetTrackStatus() <<
","
189 << sensitiveDetectorName <<
","
190 << track->GetMaterial()->GetName() <<
","
191 << step->IsFirstStepInVolume() <<
","
192 << step->IsLastStepInVolume() <<
","
193 << step->GetTotalEnergyDeposit() <<
","
194 << postStepPoint->GetProcessDefinedStep()->GetProcessType() <<
","
195 << postStepPoint->GetProcessDefinedStep()->GetProcessName() <<
","
196 << preStepPoint->GetPosition().x() <<
","
197 << preStepPoint->GetPosition().y() <<
","
198 << preStepPoint->GetPosition().z() <<
","
199 << preStepPoint->GetGlobalTime() <<
","
200 << preStepPoint->GetMomentum().x() <<
","
201 << preStepPoint->GetMomentum().y() <<
","
202 << preStepPoint->GetMomentum().z() <<
","
203 << preStepPoint->GetTotalEnergy() <<
","
204 << postStepPoint->GetPosition().x() <<
","
205 << postStepPoint->GetPosition().y() <<
","
206 << postStepPoint->GetPosition().z() <<
","
207 << postStepPoint->GetGlobalTime() <<
","
208 << postStepPoint->GetMomentum().x() <<
","
209 << postStepPoint->GetMomentum().y() <<
","
210 << postStepPoint->GetMomentum().z() <<
","
211 << postStepPoint->GetTotalEnergy() << std::endl;
bool getWriteSimSteps() const
Get the flag for writing the simulation steps into an output csv file.
static Environment & Instance()
Static method to get a reference to the Environment instance.
Class to save the full simulated trajectory of a particle.
void addPoint(float x, float y, float z, float px, float py, float pz)
Add a point to the trajectory.
bool empty() const
return true if size()==0
std::ofstream * getVREventStream() const
This method gets the output stream for the event-history steps.
The run manager controls the flow of the Geant4 program and manages the event loop(s) within a run.
UserInfo class which is used to attach additional information to Geant4 particles and tracks.
static const double mm
[millimeters]
static const double MeV
[megaelectronvolt]
static const double cm
Standard units with the value = 1.
Abstract base class for different kinds of events.