CMS 3D CMS Logo

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