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