CMS 3D CMS Logo

Phase2SteppingAction.cc
Go to the documentation of this file.
2 
5 
6 #include "G4LogicalVolumeStore.hh"
7 #include "G4ParticleTable.hh"
8 #include "G4PhysicalVolumeStore.hh"
9 #include "G4RegionStore.hh"
10 #include "G4UnitsTable.hh"
11 #include "G4SystemOfUnits.hh"
12 
15 
17  : steppingVerbose(sv), hasWatcher(hasW) {
18  theCriticalEnergyForVacuum = (p.getParameter<double>("CriticalEnergyForVacuum") * CLHEP::MeV);
19  if (0.0 < theCriticalEnergyForVacuum) {
20  killBeamPipe = true;
21  }
22  theCriticalDensity = (p.getParameter<double>("CriticalDensity") * CLHEP::g / CLHEP::cm3);
23  maxZCentralCMS = p.getParameter<double>("MaxZCentralCMS") * CLHEP::m;
24  maxTrackTime = p.getParameter<double>("MaxTrackTime") * CLHEP::ns;
25  maxTrackTimeForward = p.getParameter<double>("MaxTrackTimeForward") * CLHEP::ns;
26  maxTrackTimes = p.getParameter<std::vector<double> >("MaxTrackTimes");
27  maxTimeNames = p.getParameter<std::vector<std::string> >("MaxTimeNames");
28  deadRegionNames = p.getParameter<std::vector<std::string> >("DeadRegions");
29  maxNumberOfSteps = p.getParameter<int>("MaxNumberOfSteps");
30  ekinMins = p.getParameter<std::vector<double> >("EkinThresholds");
31  ekinNames = p.getParameter<std::vector<std::string> >("EkinNames");
32  ekinParticles = p.getParameter<std::vector<std::string> >("EkinParticles");
33 
34  edm::LogVerbatim("SimG4CoreApplication")
35  << "Phase2SteppingAction:: KillBeamPipe = " << killBeamPipe
36  << " CriticalDensity = " << theCriticalDensity * CLHEP::cm3 / CLHEP::g << " g/cm3\n"
37  << " CriticalEnergyForVacuum = " << theCriticalEnergyForVacuum / CLHEP::MeV << " Mev;"
38  << " MaxTrackTime = " << maxTrackTime / CLHEP::ns << " ns;"
39  << " MaxZCentralCMS = " << maxZCentralCMS / CLHEP::m << " m"
40  << " MaxTrackTimeForward = " << maxTrackTimeForward / CLHEP::ns << " ns"
41  << " MaxNumberOfSteps = " << maxNumberOfSteps;
42 
43  numberTimes = maxTrackTimes.size();
44  if (numberTimes > 0) {
45  for (unsigned int i = 0; i < numberTimes; i++) {
46  edm::LogVerbatim("SimG4CoreApplication")
47  << "Phase2SteppingAction::MaxTrackTime for " << maxTimeNames[i] << " is " << maxTrackTimes[i] << " ns ";
48  maxTrackTimes[i] *= ns;
49  }
50  }
51 
53  if (ndeadRegions > 0) {
54  edm::LogVerbatim("SimG4CoreApplication")
55  << "Phase2SteppingAction: Number of DeadRegions where all trackes are killed " << ndeadRegions;
56  for (unsigned int i = 0; i < ndeadRegions; ++i) {
57  edm::LogVerbatim("SimG4CoreApplication")
58  << "Phase2SteppingAction: DeadRegion " << i << ". " << deadRegionNames[i];
59  }
60  }
61  numberEkins = ekinNames.size();
62  numberPart = ekinParticles.size();
63  if (0 == numberPart) {
64  numberEkins = 0;
65  }
66 
67  if (numberEkins > 0) {
68  edm::LogVerbatim("SimG4CoreApplication")
69  << "Phase2SteppingAction::Kill following " << numberPart << " particles in " << numberEkins << " volumes";
70  for (unsigned int i = 0; i < numberPart; ++i) {
71  edm::LogVerbatim("SimG4CoreApplication") << "Phase2SteppingAction::Particle " << i << " " << ekinParticles[i]
72  << " Threshold = " << ekinMins[i] << " MeV";
73  ekinMins[i] *= CLHEP::MeV;
74  }
75  for (unsigned int i = 0; i < numberEkins; ++i) {
76  edm::LogVerbatim("SimG4CoreApplication") << "Phase2SteppingAction::LogVolume[" << i << "] = " << ekinNames[i];
77  }
78  }
79 }
80 
81 void Phase2SteppingAction::UserSteppingAction(const G4Step* aStep) {
82  if (!initialized) {
84  }
85 
86  m_g4StepSignal(aStep);
87 
88  G4Track* theTrack = aStep->GetTrack();
89  TrackStatus tstat = (theTrack->GetTrackStatus() == fAlive) ? sAlive : sKilledByProcess;
90 
91  if (theTrack->GetKineticEnergy() < 0.0) {
92  if (nWarnings < 2) {
93  ++nWarnings;
94  edm::LogWarning("SimG4CoreApplication")
95  << "Phase2SteppingAction::UserPhase2SteppingAction: Track #" << theTrack->GetTrackID() << " "
96  << theTrack->GetDefinition()->GetParticleName() << " Ekin(MeV)= " << theTrack->GetKineticEnergy() / MeV;
97  }
98  theTrack->SetKineticEnergy(0.0);
99  }
100 
101  const G4StepPoint* preStep = aStep->GetPreStepPoint();
102  const G4StepPoint* postStep = aStep->GetPostStepPoint();
103 
104  // the track is killed by the process
105  if (tstat == sKilledByProcess) {
106  if (nullptr != steppingVerbose) {
107  steppingVerbose->nextStep(aStep, fpSteppingManager, false);
108  }
109  return;
110  }
111 
112  if (sAlive == tstat && theTrack->GetCurrentStepNumber() > maxNumberOfSteps) {
113  tstat = sNumberOfSteps;
114  if (nWarnings < 5) {
115  ++nWarnings;
116  edm::LogWarning("SimG4CoreApplication")
117  << "Track #" << theTrack->GetTrackID() << " " << theTrack->GetDefinition()->GetParticleName()
118  << " E(MeV)= " << preStep->GetKineticEnergy() / MeV << " Nstep= " << theTrack->GetCurrentStepNumber()
119  << " is killed due to limit on number of steps;/n PV: " << preStep->GetPhysicalVolume()->GetName() << " at "
120  << theTrack->GetPosition() << " StepLen(mm)= " << aStep->GetStepLength();
121  }
122  }
123  const double time = theTrack->GetGlobalTime();
124 
125  // check Z-coordinate
126  if (sAlive == tstat && std::abs(theTrack->GetPosition().z()) >= maxZCentralCMS) {
128  }
129 
130  // check G4Region
131  if (sAlive == tstat) {
132  // next logical volume and next region
133  const G4LogicalVolume* lv = postStep->GetPhysicalVolume()->GetLogicalVolume();
134  const G4Region* theRegion = lv->GetRegion();
135 
136  // kill in dead regions
137  if (isInsideDeadRegion(theRegion))
138  tstat = sDeadRegion;
139 
140  // kill out of time
141  if (sAlive == tstat) {
142  if (isOutOfTimeWindow(theRegion, time))
143  tstat = sOutOfTime;
144  }
145 
146  // kill low-energy in volumes on demand
147  if (sAlive == tstat && numberEkins > 0) {
148  if (isLowEnergy(lv, theTrack))
149  tstat = sLowEnergy;
150  }
151 
152  // kill low-energy in vacuum
153  if (sAlive == tstat && killBeamPipe) {
154  if (theTrack->GetKineticEnergy() < theCriticalEnergyForVacuum &&
155  theTrack->GetDefinition()->GetPDGCharge() != 0.0 && lv->GetMaterial()->GetDensity() <= theCriticalDensity) {
156  tstat = sLowEnergyInVacuum;
157  }
158  }
159  }
160  // check transition tracker/btl and tracker/calo
161  bool isKilled = false;
162  if (sAlive == tstat || sVeryForward == tstat) {
163  // store TrackInformation about transition from one envelope to another
164  if (preStep->GetPhysicalVolume() == tracker && postStep->GetPhysicalVolume() == btl) {
165  // store transition tracker -> BTL only for tracks entering BTL for the first time
166  TrackInformation* trkinfo = static_cast<TrackInformation*>(theTrack->GetUserInformation());
167  if (!trkinfo->isFromTtoBTL() && !trkinfo->isFromBTLtoT()) {
168  trkinfo->setFromTtoBTL();
169 #ifdef EDM_ML_DEBUG
170  LogDebug("SimG4CoreApplication") << "Setting flag for Tracker -> BTL " << trkinfo->isFromTtoBTL()
171  << " IdAtBTLentrance = " << trkinfo->mcTruthID();
172 #endif
173  } else {
174  trkinfo->setBTLlooper();
175 #ifdef EDM_ML_DEBUG
176  LogDebug("SimG4CoreApplication") << "Setting flag for BTL looper " << trkinfo->isBTLlooper();
177 #endif
178  }
179  } else if (preStep->GetPhysicalVolume() == btl && postStep->GetPhysicalVolume() == tracker) {
180  // store transition BTL -> tracker
181  TrackInformation* trkinfo = static_cast<TrackInformation*>(theTrack->GetUserInformation());
182  if (!trkinfo->isFromBTLtoT()) {
183  trkinfo->setFromBTLtoT();
184 #ifdef EDM_ML_DEBUG
185  LogDebug("SimG4CoreApplication") << "Setting flag for BTL -> Tracker " << trkinfo->isFromBTLtoT();
186 #endif
187  }
188  } else if (preStep->GetPhysicalVolume() == tracker && postStep->GetPhysicalVolume() == calo) {
189  // store transition tracker -> calo
190  TrackInformation* trkinfo = static_cast<TrackInformation*>(theTrack->GetUserInformation());
191  if (!trkinfo->crossedBoundary()) {
192  trkinfo->setCrossedBoundary(theTrack);
193  }
194  } else if (preStep->GetPhysicalVolume() == calo && postStep->GetPhysicalVolume() != calo) {
195  bool backscattering(false);
196  if (postStep->GetPhysicalVolume() == tracker) {
197  backscattering = true;
198  } else if (postStep->GetPhysicalVolume() == cmse) {
199  // simple protection to avoid possible steps from calo towards the outer part of the detector, if allowed by geometry
200  // to be removed as soon as tracker-calo boundary becomes again the default
201  if (preStep->GetPosition().mag2() > postStep->GetPosition().mag2()) {
202  backscattering = true;
203  }
204  }
205  // store transition calo -> cmse to tag backscattering
206  if (backscattering) {
207  TrackInformation* trkinfo = static_cast<TrackInformation*>(theTrack->GetUserInformation());
208  if (!trkinfo->isInTrkFromBackscattering()) {
209  trkinfo->setInTrkFromBackscattering();
210 #ifdef EDM_ML_DEBUG
211  LogDebug("SimG4CoreApplication")
212  << "Setting flag for backscattering from CALO " << trkinfo->isInTrkFromBackscattering();
213 #endif
214  }
215  }
216  }
217  } else {
218  theTrack->SetTrackStatus(fStopAndKill);
219  isKilled = true;
220 #ifdef EDM_ML_DEBUG
221  PrintKilledTrack(theTrack, tstat);
222 #endif
223  }
224  if (nullptr != steppingVerbose) {
225  steppingVerbose->nextStep(aStep, fpSteppingManager, isKilled);
226  }
227 }
228 
229 bool Phase2SteppingAction::isLowEnergy(const G4LogicalVolume* lv, const G4Track* theTrack) const {
230  const double ekin = theTrack->GetKineticEnergy();
231  int pCode = theTrack->GetDefinition()->GetPDGEncoding();
232 
233  for (auto& vol : ekinVolumes) {
234  if (lv == vol) {
235  for (unsigned int i = 0; i < numberPart; ++i) {
236  if (pCode == ekinPDG[i]) {
237  return (ekin <= ekinMins[i]);
238  }
239  }
240  break;
241  }
242  }
243  return false;
244 }
245 
247  const G4PhysicalVolumeStore* pvs = G4PhysicalVolumeStore::GetInstance();
248  for (auto const& pvcite : *pvs) {
249  const G4String& pvname = pvcite->GetName();
250  if (pvname == "Tracker" || pvname == "tracker:Tracker_1") {
251  tracker = pvcite;
252  } else if (pvname == "CALO" || pvname == "caloBase:CALO_1") {
253  calo = pvcite;
254  } else if (pvname == "BarrelTimingLayer" || pvname == "btl:BarrelTimingLayer_1") {
255  btl = pvcite;
256  } else if (pvname == "CMSE" || pvname == "cms:CMSE_1") {
257  cmse = pvcite;
258  }
259  if (tracker && calo && btl && cmse)
260  break;
261  }
262  edm::LogVerbatim("SimG4CoreApplication")
263  << "Phase2SteppingAction: pointer for Tracker " << tracker << " and for Calo " << calo << " and for BTL " << btl;
264 
265  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
266  if (numberEkins > 0) {
267  ekinVolumes.resize(numberEkins, nullptr);
268  for (auto const& lvcite : *lvs) {
269  const G4String& lvname = lvcite->GetName();
270  for (unsigned int i = 0; i < numberEkins; ++i) {
271  if (lvname == (G4String)(ekinNames[i])) {
272  ekinVolumes[i] = lvcite;
273  break;
274  }
275  }
276  }
277  for (unsigned int i = 0; i < numberEkins; ++i) {
278  edm::LogVerbatim("SimG4CoreApplication") << ekinVolumes[i]->GetName() << " with pointer " << ekinVolumes[i];
279  }
280  }
281 
282  if (numberPart > 0) {
283  G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
284  ekinPDG.resize(numberPart, 0);
285  for (unsigned int i = 0; i < numberPart; ++i) {
286  const G4ParticleDefinition* part = theParticleTable->FindParticle(ekinParticles[i]);
287  if (nullptr != part)
288  ekinPDG[i] = part->GetPDGEncoding();
289  edm::LogVerbatim("SimG4CoreApplication") << "Particle " << ekinParticles[i] << " with PDG code " << ekinPDG[i]
290  << " and KE cut off " << ekinMins[i] / MeV << " MeV";
291  }
292  }
293 
294  const G4RegionStore* rs = G4RegionStore::GetInstance();
295  if (numberTimes > 0) {
296  maxTimeRegions.resize(numberTimes, nullptr);
297  for (auto const& rcite : *rs) {
298  const G4String& rname = rcite->GetName();
299  for (unsigned int i = 0; i < numberTimes; ++i) {
300  if (rname == (G4String)(maxTimeNames[i])) {
301  maxTimeRegions[i] = rcite;
302  break;
303  }
304  }
305  }
306  }
307  if (ndeadRegions > 0) {
308  deadRegions.resize(ndeadRegions, nullptr);
309  for (auto const& rcite : *rs) {
310  const G4String& rname = rcite->GetName();
311  for (unsigned int i = 0; i < ndeadRegions; ++i) {
312  if (rname == (G4String)(deadRegionNames[i])) {
313  deadRegions[i] = rcite;
314  break;
315  }
316  }
317  }
318  }
319  return true;
320 }
321 
322 void Phase2SteppingAction::PrintKilledTrack(const G4Track* aTrack, const TrackStatus& tst) const {
323  std::string vname = "";
324  std::string rname = "";
325  std::string typ = " ";
326  switch (tst) {
327  case sDeadRegion:
328  typ = " in dead region ";
329  break;
330  case sOutOfTime:
331  typ = " out of time window ";
332  break;
333  case sLowEnergy:
334  typ = " low energy limit ";
335  break;
336  case sLowEnergyInVacuum:
337  typ = " low energy limit in vacuum ";
338  break;
339  case sEnergyDepNaN:
340  typ = " energy deposition is NaN ";
341  break;
342  case sVeryForward:
343  typ = " very forward track ";
344  break;
345  case sNumberOfSteps:
346  typ = " too many steps ";
347  break;
348  default:
349  break;
350  }
351  G4VPhysicalVolume* pv = aTrack->GetNextVolume();
352  vname = pv->GetLogicalVolume()->GetName();
353  rname = pv->GetLogicalVolume()->GetRegion()->GetName();
354 
355  const double ekin = aTrack->GetKineticEnergy();
356  if (ekin < 2 * CLHEP::MeV) {
357  return;
358  }
359  edm::LogWarning("SimG4CoreApplication")
360  << "Track #" << aTrack->GetTrackID() << " StepN= " << aTrack->GetCurrentStepNumber() << " "
361  << aTrack->GetDefinition()->GetParticleName() << " E(MeV)=" << ekin / CLHEP::MeV
362  << " T(ns)=" << aTrack->GetGlobalTime() / CLHEP::ns << " is killed due to " << typ << "\n LV: " << vname << " ("
363  << rname << ") at " << aTrack->GetPosition() << " step(cm)=" << aTrack->GetStep()->GetStepLength() / CLHEP::cm;
364 }
Log< level::Info, true > LogVerbatim
bool crossedBoundary() const
const CMSSteppingVerbose * steppingVerbose
std::vector< std::string > maxTimeNames
bool isOutOfTimeWindow(const G4Region *reg, const double &time) const
void setInTrkFromBackscattering()
SimActivityRegistry::G4StepSignal m_g4StepSignal
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
const G4VPhysicalVolume * calo
Phase2SteppingAction(const CMSSteppingVerbose *, const edm::ParameterSet &, bool hasW)
std::vector< std::string > ekinNames
const G4VPhysicalVolume * tracker
std::vector< int > ekinPDG
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isInTrkFromBackscattering() const
std::vector< G4LogicalVolume * > ekinVolumes
const G4VPhysicalVolume * btl
bool isFromBTLtoT() const
std::vector< std::string > ekinParticles
int mcTruthID() const
part
Definition: HCALResponse.h:20
void PrintKilledTrack(const G4Track *, const TrackStatus &) const
bool isBTLlooper() const
std::vector< double > ekinMins
const G4String rname[NREG]
std::vector< std::string > deadRegionNames
std::vector< const G4Region * > maxTimeRegions
bool isLowEnergy(const G4LogicalVolume *, const G4Track *) const
bool isInsideDeadRegion(const G4Region *reg) const
void nextStep(const G4Step *, const G4SteppingManager *ptr, bool isKilled) const
bool isFromTtoBTL() const
std::vector< const G4Region * > deadRegions
Log< level::Warning, false > LogWarning
void setCrossedBoundary(const G4Track *track)
TrackStatus
Definition: Common.h:9
const G4VPhysicalVolume * cmse
std::vector< double > maxTrackTimes
void UserSteppingAction(const G4Step *aStep) final
#define LogDebug(id)