Belle II Software development
SteppingAction Class Reference

The Class for the stepping action. More...

#include <SteppingAction.h>

Inheritance diagram for SteppingAction:

Public Member Functions

 SteppingAction ()
 Constructor.
 
virtual ~SteppingAction ()
 Destructor.
 
void setMaxNumberSteps (int maxSteps)
 Sets the maximum number of steps before a track is stopped and killed.
 
void setStoreTrajectories (bool store)
 Sets the trajectory option to enable storing of the simulated particle trajectories.
 
void setAbsorbersR (const std::vector< float > &vec)
 Sets the radii of absorbers for killing tracks across them.
 
virtual void UserSteppingAction (const G4Step *step)
 The method will be called at each step during simulation.
 

Protected Attributes

int m_maxNumberSteps
 The maximum number of steps before the track transportation is stopped and the track is killed.
 
bool m_storeTrajectories = false
 if true, check if the track has attached trajectory info and append step information if necessary
 
std::vector< float > m_absorbers
 The absorbers defined at given radii where tracks across them will be destroyed.
 

Private Member Functions

void writeVREventStep (const G4Step *, const G4Track *)
 Method to write (almost) each G4Step to the VR event file.
 

Private Attributes

bool m_writeSimSteps {false}
 Flag for writing out the simulation steps.
 

Detailed Description

The Class for the stepping action.

This class is called for each step during the Geant4 transportation. It implements protection mechanisms to remove unreasonable tracks.

Definition at line 31 of file SteppingAction.h.

Constructor & Destructor Documentation

◆ SteppingAction()

Constructor.

Definition at line 27 of file SteppingAction.cc.

28{
29 //Default value for the maximum number of steps
30 m_maxNumberSteps = 100000;
31 if (false) {
32 G4Step* aStep;
33 UserSteppingAction(aStep);
34 }
36}
bool getWriteSimSteps() const
Get the flag for writing the simulation steps into an output csv file.
Definition: Environment.h:219
static Environment & Instance()
Static method to get a reference to the Environment instance.
Definition: Environment.cc:28
bool m_writeSimSteps
Flag for writing out the simulation steps.
virtual void UserSteppingAction(const G4Step *step)
The method will be called at each step during simulation.
int m_maxNumberSteps
The maximum number of steps before the track transportation is stopped and the track is killed.

◆ ~SteppingAction()

~SteppingAction ( )
virtual

Destructor.

Definition at line 39 of file SteppingAction.cc.

40{
41
42}

Member Function Documentation

◆ setAbsorbersR()

void setAbsorbersR ( const std::vector< float > &  vec)
inline

Sets the radii of absorbers for killing tracks across them.

This set is used in the PXD only simulation for PXD gain calibration.

Parameters
vecThe c++ vector with the radii of absorbers in cm

Definition at line 57 of file SteppingAction.h.

57{ m_absorbers = vec; };
std::vector< float > m_absorbers
The absorbers defined at given radii where tracks across them will be destroyed.

◆ setMaxNumberSteps()

void setMaxNumberSteps ( int  maxSteps)
inline

Sets the maximum number of steps before a track is stopped and killed.

Parameters
maxStepsThe maximum number of steps.

Definition at line 48 of file SteppingAction.h.

48{ m_maxNumberSteps = maxSteps; };

◆ setStoreTrajectories()

void setStoreTrajectories ( bool  store)
inline

Sets the trajectory option to enable storing of the simulated particle trajectories.

Definition at line 51 of file SteppingAction.h.

51{ m_storeTrajectories = store; }
bool m_storeTrajectories
if true, check if the track has attached trajectory info and append step information if necessary

◆ UserSteppingAction()

void UserSteppingAction ( const G4Step *  step)
virtual

The method will be called at each step during simulation.

Parameters
stepThe pointer of current step.

Definition at line 45 of file SteppingAction.cc.

46{
47 G4Track* track = step->GetTrack();
48
49 //------------------------------
50 // Check for NULL world volume
51 //------------------------------
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);
57 return;
58 }
59
60 //------------------------------
61 // Check for absorbers
62 //------------------------------
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)) {
67 //B2WARNING("SteppingAction: Track across absorbers, terminating!\n"
68 //<< "step_no=" << track->GetCurrentStepNumber() << " type=" << track->GetDefinition()->GetParticleName()
69 //<< "\n position=" << G4BestUnit(track->GetPosition(), "Length") << " momentum=" << G4BestUnit(track->GetMomentum(), "Energy") << "\n PrePos.perp=" << stepPrePos.perp() << ", PostPos.perp=" << stepPostPos.perp() << " cm" );
70 track->SetTrackStatus(fStopAndKill);
71 return;
72 }
73 }
74
75 // If we are running this job for producing virtual reality events, let's run the relevant method
76 if (m_writeSimSteps) {
77 writeVREventStep(step, track);
78 }
79
80 //---------------------------------------
81 // Check for very high number of steps.
82 //---------------------------------------
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);
88 return;
89 }
90
91 //-----------------------------------------------------------
92 // Check if there is an attached trajectory. If so, fill it.
93 //-----------------------------------------------------------
95 TrackInfo* info = dynamic_cast<TrackInfo*>(track->GetUserInformation());
96 if (info && info->getTrajectory()) {
97 MCParticleTrajectory& trajectory = *(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;
101 trajectory.addPoint(
102 stepPos.x(), stepPos.y(), stepPos.z(),
103 stepMom.x(), stepMom.y(), stepMom.z()
104 );
105 }
106 const G4ThreeVector stepPos = step->GetPostStepPoint()->GetPosition() / CLHEP::mm * Unit::mm;
107 const G4ThreeVector stepMom = step->GetPostStepPoint()->GetMomentum() / CLHEP::MeV * Unit::MeV;
108 trajectory.addPoint(
109 stepPos.x(), stepPos.y(), stepPos.z(),
110 stepMom.x(), stepMom.y(), stepMom.z()
111 );
112 }
113 }
114}
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
void writeVREventStep(const G4Step *, const G4Track *)
Method to write (almost) each G4Step to the VR event file.
static const double mm
[millimeters]
Definition: Unit.h:70
static const double MeV
[megaelectronvolt]
Definition: Unit.h:114
static const double cm
Standard units with the value = 1.
Definition: Unit.h:47

◆ writeVREventStep()

void writeVREventStep ( const G4Step *  step,
const G4Track *  track 
)
private

Method to write (almost) each G4Step to the VR event file.

Definition at line 117 of file SteppingAction.cc.

118{
119 // There must be an open output file (opened in simulation's EventAction)
120 RunManager& runManager = RunManager::Instance();
121 const EventAction* eventAction = (const Belle2::Simulation::EventAction*)runManager.GetUserEventAction();
122 if (eventAction == nullptr)
123 return;
124 std::ofstream* output = eventAction->getVREventStream();
125 if (output == nullptr)
126 return;
127 if (!output->is_open())
128 return;
129 // Ignore (hyper)nuclei heavier than alphas
130 if (abs(track->GetDefinition()->GetPDGEncoding()) > 1000020040)
131 return;
132 // Limit the VR event history to 100 ns
133 G4StepPoint* postStepPoint = step->GetPostStepPoint();
134 if (postStepPoint->GetGlobalTime() > 100.0)
135 return;
136 // Discard soft particles (KE < 500 keV) unless they are optical photons (for which KE=0)
137 G4StepPoint* preStepPoint = step->GetPreStepPoint();
138 double KE = preStepPoint->GetTotalEnergy() - track->GetDefinition()->GetPDGMass();
139 if ((KE < 0.0005) && (track->GetDefinition()->GetParticleName() != "opticalphoton"))
140 return;
141
142 // We will definitely write one record to the VR event file
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";
170 }
171 } else if (pVolName.compare(0, 15, "StripSensitive_") == 0) {
172 sensitiveDetectorName = "EKLM";
173 }
174// Content of each record:
175// TrackID,ParentID,ParticleName,Mass,Charge,StepNumber,Status,VolumeName,
176// MaterialName,IsFirstStepInVolume,IsLastStepInVolume,EnergyDeposit,
177// ProcessType,ProcessName,PrePointX,PrePointY,PrePointZ,PrePointT,
178// PrePointPX,PrePointPY,PrePointPZ,PrePointE,PostPointX,PostPointY,
179// PostPointZ,PostPointT,PostPointPX,PostPointPY,PostPointPZ,PostPointE
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() << ","
188 << pVolName << ","
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;
212}
The Event Action class.
Definition: EventAction.h:36
static RunManager & Instance()
Static method to get a reference to the RunManager instance.
Definition: RunManager.cc:29

Member Data Documentation

◆ m_absorbers

std::vector<float> m_absorbers
protected

The absorbers defined at given radii where tracks across them will be destroyed.

Definition at line 70 of file SteppingAction.h.

◆ m_maxNumberSteps

int m_maxNumberSteps
protected

The maximum number of steps before the track transportation is stopped and the track is killed.

Definition at line 67 of file SteppingAction.h.

◆ m_storeTrajectories

bool m_storeTrajectories = false
protected

if true, check if the track has attached trajectory info and append step information if necessary

Definition at line 69 of file SteppingAction.h.

◆ m_writeSimSteps

bool m_writeSimSteps {false}
private

Flag for writing out the simulation steps.

Definition at line 78 of file SteppingAction.h.


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