CMS 3D CMS Logo

SteppingAction.cc
Go to the documentation of this file.
4 
5 #include "G4LogicalVolumeStore.hh"
6 #include "G4ParticleTable.hh"
7 #include "G4PhysicalVolumeStore.hh"
8 #include "G4RegionStore.hh"
9 #include "G4UnitsTable.hh"
10 #include "G4SystemOfUnits.hh"
11 
14 
15 //#define DebugLog
16 
18  : eventAction_(e),
19  tracker(nullptr),
20  calo(nullptr),
21  steppingVerbose(sv),
22  nWarnings(0),
23  initialized(false),
24  killBeamPipe(false),
25  hasWatcher(hasW) {
26  theCriticalEnergyForVacuum = (p.getParameter<double>("CriticalEnergyForVacuum") * CLHEP::MeV);
27  if (0.0 < theCriticalEnergyForVacuum) {
28  killBeamPipe = true;
29  }
30  m_CMStoZDCtransport = (p.getParameter<bool>("CMStoZDCtransport"));
31  theCriticalDensity = (p.getParameter<double>("CriticalDensity") * CLHEP::g / CLHEP::cm3);
32  maxZCentralCMS = p.getParameter<double>("MaxZCentralCMS") * CLHEP::m;
33  maxTrackTime = p.getParameter<double>("MaxTrackTime") * CLHEP::ns;
34  maxTrackTimeForward = p.getParameter<double>("MaxTrackTimeForward") * CLHEP::ns;
35  maxTrackTimes = p.getParameter<std::vector<double> >("MaxTrackTimes");
36  maxTimeNames = p.getParameter<std::vector<std::string> >("MaxTimeNames");
37  deadRegionNames = p.getParameter<std::vector<std::string> >("DeadRegions");
38  maxNumberOfSteps = p.getParameter<int>("MaxNumberOfSteps");
39  ekinMins = p.getParameter<std::vector<double> >("EkinThresholds");
40  ekinNames = p.getParameter<std::vector<std::string> >("EkinNames");
41  ekinParticles = p.getParameter<std::vector<std::string> >("EkinParticles");
42 
43  edm::LogVerbatim("SimG4CoreApplication")
44  << "SteppingAction:: KillBeamPipe = " << killBeamPipe
45  << " CriticalDensity = " << theCriticalDensity * CLHEP::cm3 / CLHEP::g << " g/cm3\n"
46  << " CriticalEnergyForVacuum = " << theCriticalEnergyForVacuum / CLHEP::MeV << " Mev;"
47  << " MaxTrackTime = " << maxTrackTime / CLHEP::ns << " ns;"
48  << " MaxZCentralCMS = " << maxZCentralCMS / CLHEP::m << " m"
49  << " MaxTrackTimeForward = " << maxTrackTimeForward / CLHEP::ns << " ns"
50  << " MaxNumberOfSteps = " << maxNumberOfSteps << " ZDC: " << m_CMStoZDCtransport;
51 
52  numberTimes = maxTrackTimes.size();
53  if (numberTimes > 0) {
54  for (unsigned int i = 0; i < numberTimes; i++) {
55  edm::LogVerbatim("SimG4CoreApplication")
56  << "SteppingAction::MaxTrackTime for " << maxTimeNames[i] << " is " << maxTrackTimes[i] << " ns ";
57  maxTrackTimes[i] *= ns;
58  }
59  }
60 
62  if (ndeadRegions > 0) {
63  edm::LogVerbatim("SimG4CoreApplication")
64  << "SteppingAction: Number of DeadRegions where all trackes are killed " << ndeadRegions;
65  for (unsigned int i = 0; i < ndeadRegions; ++i) {
66  edm::LogVerbatim("SimG4CoreApplication") << "SteppingAction: DeadRegion " << i << ". " << deadRegionNames[i];
67  }
68  }
69  numberEkins = ekinNames.size();
70  numberPart = ekinParticles.size();
71  if (0 == numberPart) {
72  numberEkins = 0;
73  }
74 
75  if (numberEkins > 0) {
76  edm::LogVerbatim("SimG4CoreApplication")
77  << "SteppingAction::Kill following " << numberPart << " particles in " << numberEkins << " volumes";
78  for (unsigned int i = 0; i < numberPart; ++i) {
79  edm::LogVerbatim("SimG4CoreApplication")
80  << "SteppingAction::Particle " << i << " " << ekinParticles[i] << " Threshold = " << ekinMins[i] << " MeV";
81  ekinMins[i] *= CLHEP::MeV;
82  }
83  for (unsigned int i = 0; i < numberEkins; ++i) {
84  edm::LogVerbatim("SimG4CoreApplication") << "SteppingAction::LogVolume[" << i << "] = " << ekinNames[i];
85  }
86  }
87 }
88 
90 
91 void SteppingAction::UserSteppingAction(const G4Step* aStep) {
92  if (!initialized) {
94  }
95 
96  //if(hasWatcher) { m_g4StepSignal(aStep); }
97  m_g4StepSignal(aStep);
98 
99  G4Track* theTrack = aStep->GetTrack();
100  TrackStatus tstat = (theTrack->GetTrackStatus() == fAlive) ? sAlive : sKilledByProcess;
101 
102  if (theTrack->GetKineticEnergy() < 0.0) {
103  if (nWarnings < 2) {
104  ++nWarnings;
105  edm::LogWarning("SimG4CoreApplication")
106  << "SteppingAction::UserSteppingAction: Track #" << theTrack->GetTrackID() << " "
107  << theTrack->GetDefinition()->GetParticleName() << " Ekin(MeV)= " << theTrack->GetKineticEnergy() / MeV;
108  }
109  theTrack->SetKineticEnergy(0.0);
110  }
111 
112  const G4StepPoint* preStep = aStep->GetPreStepPoint();
113  const G4StepPoint* postStep = aStep->GetPostStepPoint();
114 
115  // NaN energy deposit
116  if (edm::isNotFinite(aStep->GetTotalEnergyDeposit())) {
117  tstat = sEnergyDepNaN;
118  if (nWarnings < 5) {
119  ++nWarnings;
120  edm::LogWarning("SimG4CoreApplication")
121  << "Track #" << theTrack->GetTrackID() << " " << theTrack->GetDefinition()->GetParticleName()
122  << " E(MeV)= " << preStep->GetKineticEnergy() / MeV << " Nstep= " << theTrack->GetCurrentStepNumber()
123  << " is killed due to edep=NaN inside PV: " << preStep->GetPhysicalVolume()->GetName() << " at "
124  << theTrack->GetPosition() << " StepLen(mm)= " << aStep->GetStepLength();
125  }
126  }
127 
128  // the track is killed by the process
129  if (tstat == sKilledByProcess) {
130  if (nullptr != steppingVerbose) {
131  steppingVerbose->NextStep(aStep, fpSteppingManager, false);
132  }
133  return;
134  }
135 
136  if (sAlive == tstat && theTrack->GetCurrentStepNumber() > maxNumberOfSteps) {
137  tstat = sNumberOfSteps;
138  if (nWarnings < 5) {
139  ++nWarnings;
140  edm::LogWarning("SimG4CoreApplication")
141  << "Track #" << theTrack->GetTrackID() << " " << theTrack->GetDefinition()->GetParticleName()
142  << " E(MeV)= " << preStep->GetKineticEnergy() / MeV << " Nstep= " << theTrack->GetCurrentStepNumber()
143  << " is killed due to limit on number of steps;/n PV: " << preStep->GetPhysicalVolume()->GetName() << " at "
144  << theTrack->GetPosition() << " StepLen(mm)= " << aStep->GetStepLength();
145  }
146  }
147  const double time = theTrack->GetGlobalTime();
148 
149  // check Z-coordinate
150  if (sAlive == tstat && std::abs(theTrack->GetPosition().z()) >= maxZCentralCMS) {
152  }
153 
154  // check G4Region
155  if (sAlive == tstat) {
156  // next logical volume and next region
157  const G4LogicalVolume* lv = postStep->GetPhysicalVolume()->GetLogicalVolume();
158  const G4Region* theRegion = lv->GetRegion();
159 
160  // kill in dead regions
161  if (isInsideDeadRegion(theRegion) && !isForZDC(lv, std::abs(theTrack->GetParticleDefinition()->GetPDGEncoding()))) {
162  tstat = sDeadRegion;
163  }
164 
165  // kill out of time
166  if (sAlive == tstat) {
167  if (isOutOfTimeWindow(theRegion, time))
168  tstat = sOutOfTime;
169  }
170 
171  // kill low-energy in volumes on demand
172  if (sAlive == tstat && numberEkins > 0) {
173  if (isLowEnergy(lv, theTrack))
174  tstat = sLowEnergy;
175  }
176 
177  // kill low-energy in vacuum
178  if (sAlive == tstat && killBeamPipe) {
179  if (theTrack->GetKineticEnergy() < theCriticalEnergyForVacuum &&
180  theTrack->GetDefinition()->GetPDGCharge() != 0.0 && lv->GetMaterial()->GetDensity() <= theCriticalDensity) {
181  tstat = sLowEnergyInVacuum;
182  }
183  }
184  }
185  // check transition tracker/calo
186  bool isKilled = false;
187  if (sAlive == tstat || sVeryForward == tstat) {
188  if (preStep->GetPhysicalVolume() == tracker && postStep->GetPhysicalVolume() == calo) {
189  math::XYZVectorD pos((postStep->GetPosition()).x(), (postStep->GetPosition()).y(), (postStep->GetPosition()).z());
190 
191  math::XYZTLorentzVectorD mom((postStep->GetMomentum()).x(),
192  (postStep->GetMomentum()).y(),
193  (postStep->GetMomentum()).z(),
194  postStep->GetTotalEnergy());
195 
196  uint32_t id = theTrack->GetTrackID();
197 
198  std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> p(pos, mom);
200  }
201  } else {
202  theTrack->SetTrackStatus(fStopAndKill);
203  isKilled = true;
204 #ifdef DebugLog
205  PrintKilledTrack(theTrack, tstat);
206 #endif
207  }
208  if (nullptr != steppingVerbose) {
209  steppingVerbose->NextStep(aStep, fpSteppingManager, isKilled);
210  }
211 }
212 
213 bool SteppingAction::isLowEnergy(const G4LogicalVolume* lv, const G4Track* theTrack) const {
214  const double ekin = theTrack->GetKineticEnergy();
215  int pCode = theTrack->GetDefinition()->GetPDGEncoding();
216 
217  for (auto& vol : ekinVolumes) {
218  if (lv == vol) {
219  for (unsigned int i = 0; i < numberPart; ++i) {
220  if (pCode == ekinPDG[i]) {
221  return (ekin <= ekinMins[i]);
222  }
223  }
224  break;
225  }
226  }
227  return false;
228 }
229 
231  const G4PhysicalVolumeStore* pvs = G4PhysicalVolumeStore::GetInstance();
232  for (auto const& pvcite : *pvs) {
233  const G4String& pvname = pvcite->GetName();
234  if (pvname == "Tracker" || pvname == "tracker:Tracker_1")
235  tracker = pvcite;
236  else if (pvname == "CALO" || pvname == "caloBase:CALO_1")
237  calo = pvcite;
238 
239  if (tracker && calo)
240  break;
241  }
242  edm::LogVerbatim("SimG4CoreApplication")
243  << "SteppingAction: pointer for Tracker " << tracker << " and for Calo " << calo;
244 
245  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
246  if (numberEkins > 0) {
247  ekinVolumes.resize(numberEkins, nullptr);
248  for (auto const& lvcite : *lvs) {
249  const G4String& lvname = lvcite->GetName();
250  if (lvname == "CMStoZDC" || lvname == "cmsextent:CMStoZDC") {
251  m_CMStoZDC = lvcite;
252  }
253  for (unsigned int i = 0; i < numberEkins; ++i) {
254  if (lvname == (G4String)(ekinNames[i])) {
255  ekinVolumes[i] = lvcite;
256  break;
257  }
258  }
259  }
260  for (unsigned int i = 0; i < numberEkins; ++i) {
261  edm::LogVerbatim("SimG4CoreApplication") << ekinVolumes[i]->GetName() << " with pointer " << ekinVolumes[i];
262  }
263  }
264 
265  if (numberPart > 0) {
266  G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
267  ekinPDG.resize(numberPart, 0);
268  for (unsigned int i = 0; i < numberPart; ++i) {
269  const G4ParticleDefinition* part = theParticleTable->FindParticle(ekinParticles[i]);
270  if (nullptr != part)
271  ekinPDG[i] = part->GetPDGEncoding();
272  edm::LogVerbatim("SimG4CoreApplication") << "Particle " << ekinParticles[i] << " with PDG code " << ekinPDG[i]
273  << " and KE cut off " << ekinMins[i] / MeV << " MeV";
274  }
275  }
276 
277  const G4RegionStore* rs = G4RegionStore::GetInstance();
278  if (numberTimes > 0) {
279  maxTimeRegions.resize(numberTimes, nullptr);
280  for (auto const& rcite : *rs) {
281  const G4String& rname = rcite->GetName();
282  for (unsigned int i = 0; i < numberTimes; ++i) {
283  if (rname == (G4String)(maxTimeNames[i])) {
284  maxTimeRegions[i] = rcite;
285  break;
286  }
287  }
288  }
289  }
290  if (ndeadRegions > 0) {
291  deadRegions.resize(ndeadRegions, nullptr);
292  for (auto const& rcite : *rs) {
293  const G4String& rname = rcite->GetName();
294  for (unsigned int i = 0; i < ndeadRegions; ++i) {
295  if (rname == (G4String)(deadRegionNames[i])) {
296  deadRegions[i] = rcite;
297  break;
298  }
299  }
300  }
301  }
302  return true;
303 }
304 
305 void SteppingAction::PrintKilledTrack(const G4Track* aTrack, const TrackStatus& tst) const {
306  std::string vname = "";
307  std::string rname = "";
308  std::string typ = " ";
309  switch (tst) {
310  case sDeadRegion:
311  typ = " in dead region ";
312  break;
313  case sOutOfTime:
314  typ = " out of time window ";
315  break;
316  case sLowEnergy:
317  typ = " low energy limit ";
318  break;
319  case sLowEnergyInVacuum:
320  typ = " low energy limit in vacuum ";
321  break;
322  case sEnergyDepNaN:
323  typ = " energy deposition is NaN ";
324  break;
325  case sVeryForward:
326  typ = " very forward track ";
327  break;
328  case sNumberOfSteps:
329  typ = " too many steps ";
330  break;
331  default:
332  break;
333  }
334  G4VPhysicalVolume* pv = aTrack->GetNextVolume();
335  vname = pv->GetLogicalVolume()->GetName();
336  rname = pv->GetLogicalVolume()->GetRegion()->GetName();
337 
338  const double ekin = aTrack->GetKineticEnergy();
339  if (ekin < 2 * CLHEP::MeV) {
340  return;
341  }
342  edm::LogWarning("SimG4CoreApplication")
343  << "Track #" << aTrack->GetTrackID() << " StepN= " << aTrack->GetCurrentStepNumber() << " "
344  << aTrack->GetDefinition()->GetParticleName() << " E(MeV)=" << ekin / CLHEP::MeV
345  << " T(ns)=" << aTrack->GetGlobalTime() / CLHEP::ns << " is killed due to " << typ << "\n LV: " << vname << " ("
346  << rname << ") at " << aTrack->GetPosition() << " step(cm)=" << aTrack->GetStep()->GetStepLength() / CLHEP::cm;
347 }
bool m_CMStoZDCtransport
Log< level::Info, true > LogVerbatim
const G4LogicalVolume * m_CMStoZDC
void UserSteppingAction(const G4Step *aStep) final
const G4VPhysicalVolume * tracker
std::vector< int > ekinPDG
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
std::vector< const G4Region * > deadRegions
double theCriticalDensity
SimActivityRegistry::G4StepSignal m_g4StepSignal
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
~SteppingAction() override
unsigned int numberPart
bool isOutOfTimeWindow(const G4Region *reg, const double &time) const
const CMSSteppingVerbose * steppingVerbose
unsigned int numberEkins
bool isForZDC(const G4LogicalVolume *lv, int pdg) const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
TrackStatus
G4int maxNumberOfSteps
void PrintKilledTrack(const G4Track *, const TrackStatus &) const
std::vector< std::string > deadRegionNames
const G4VPhysicalVolume * calo
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isLowEnergy(const G4LogicalVolume *, const G4Track *) const
unsigned int numberTimes
std::vector< double > maxTrackTimes
bool isInsideDeadRegion(const G4Region *reg) const
std::vector< G4LogicalVolume * > ekinVolumes
double theCriticalEnergyForVacuum
EventAction * eventAction_
double maxZCentralCMS
std::vector< double > ekinMins
std::vector< std::string > maxTimeNames
part
Definition: HCALResponse.h:20
std::vector< std::string > ekinNames
double maxTrackTimeForward
const G4String rname[NREG]
void addTkCaloStateInfo(uint32_t t, const std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > &p)
Definition: EventAction.cc:71
void NextStep(const G4Step *, const G4SteppingManager *ptr, bool isKilled) const
SteppingAction(EventAction *ea, const edm::ParameterSet &ps, const CMSSteppingVerbose *, bool hasW)
Log< level::Warning, false > LogWarning
unsigned int ndeadRegions
std::vector< std::string > ekinParticles
Definition: Common.h:9
std::vector< const G4Region * > maxTimeRegions
unsigned int nWarnings