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 "G4SystemOfUnits.hh"
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] *= 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() << " Ekin(MeV)= " << theTrack->GetKineticEnergy() / MeV;
107  }
108  theTrack->SetKineticEnergy(0.0);
109  }
110 
111  const G4StepPoint* preStep = aStep->GetPreStepPoint();
112  const G4StepPoint* postStep = aStep->GetPostStepPoint();
113 
114  // the track is killed by the process
115  if (tstat == sKilledByProcess) {
116  if (nullptr != steppingVerbose) {
117  steppingVerbose->nextStep(aStep, fpSteppingManager, false);
118  }
119  return;
120  }
121 
122  if (sAlive == tstat && theTrack->GetCurrentStepNumber() > maxNumberOfSteps) {
123  tstat = sNumberOfSteps;
124  if (nWarnings < 5) {
125  ++nWarnings;
126  edm::LogWarning("SimG4CoreApplication")
127  << "Track #" << theTrack->GetTrackID() << " " << theTrack->GetDefinition()->GetParticleName()
128  << " E(MeV)= " << preStep->GetKineticEnergy() / MeV << " Nstep= " << theTrack->GetCurrentStepNumber()
129  << " is killed due to limit on number of steps;/n PV: " << preStep->GetPhysicalVolume()->GetName() << " at "
130  << theTrack->GetPosition() << " StepLen(mm)= " << aStep->GetStepLength();
131  }
132  }
133  const double time = theTrack->GetGlobalTime();
134 
135  // check Z-coordinate
136  if (sAlive == tstat && std::abs(theTrack->GetPosition().z()) >= maxZCentralCMS) {
138  }
139 
140  // check G4Region
141  if (sAlive == tstat) {
142  // next logical volume and next region
143  const G4LogicalVolume* lv = postStep->GetPhysicalVolume()->GetLogicalVolume();
144  const G4Region* theRegion = lv->GetRegion();
145 
146  // kill in dead regions
147  if (isInsideDeadRegion(theRegion))
148  tstat = sDeadRegion;
149 
150  // kill out of time
151  if (sAlive == tstat) {
152  if (isOutOfTimeWindow(theRegion, time))
153  tstat = sOutOfTime;
154  }
155 
156  // kill low-energy in volumes on demand
157  if (sAlive == tstat && numberEkins > 0) {
158  if (isLowEnergy(lv, theTrack))
159  tstat = sLowEnergy;
160  }
161 
162  // kill low-energy in vacuum
163  if (sAlive == tstat && killBeamPipe) {
164  if (theTrack->GetKineticEnergy() < theCriticalEnergyForVacuum &&
165  theTrack->GetDefinition()->GetPDGCharge() != 0.0 && lv->GetMaterial()->GetDensity() <= theCriticalDensity) {
166  tstat = sLowEnergyInVacuum;
167  }
168  }
169  }
170  // check transition tracker/btl and tracker/calo
171  bool isKilled = false;
172  if (sAlive == tstat || sVeryForward == tstat) {
173  // store TrackInformation about transition from one envelope to another
174  if (preStep->GetPhysicalVolume() == tracker && postStep->GetPhysicalVolume() == btl) {
175  // store transition tracker -> BTL only for tracks entering BTL for the first time
176  TrackInformation* trkinfo = static_cast<TrackInformation*>(theTrack->GetUserInformation());
177  if (!trkinfo->isFromTtoBTL() && !trkinfo->isFromBTLtoT()) {
178  trkinfo->setFromTtoBTL();
179 #ifdef EDM_ML_DEBUG
180  LogDebug("SimG4CoreApplication") << "Setting flag for Tracker -> BTL " << trkinfo->isFromTtoBTL()
181  << " IdAtBTLentrance = " << trkinfo->mcTruthID();
182 #endif
183  } else {
184  trkinfo->setBTLlooper();
185 #ifdef EDM_ML_DEBUG
186  LogDebug("SimG4CoreApplication") << "Setting flag for BTL looper " << trkinfo->isBTLlooper();
187 #endif
188  }
189  } else if (preStep->GetPhysicalVolume() == btl && postStep->GetPhysicalVolume() == tracker) {
190  // store transition BTL -> tracker
191  TrackInformation* trkinfo = static_cast<TrackInformation*>(theTrack->GetUserInformation());
192  if (!trkinfo->isFromBTLtoT()) {
193  trkinfo->setFromBTLtoT();
194 #ifdef EDM_ML_DEBUG
195  LogDebug("SimG4CoreApplication") << "Setting flag for BTL -> Tracker " << trkinfo->isFromBTLtoT();
196 #endif
197  }
198  } else if (preStep->GetPhysicalVolume() == tracker && postStep->GetPhysicalVolume() == calo) {
199  // store transition tracker -> calo
200  TrackInformation* trkinfo = static_cast<TrackInformation*>(theTrack->GetUserInformation());
201  if (!trkinfo->crossedBoundary()) {
202  trkinfo->setCrossedBoundary(theTrack);
203  }
204  } else if (preStep->GetPhysicalVolume() == calo && postStep->GetPhysicalVolume() != calo) {
205  bool backscattering(false);
206  if (postStep->GetPhysicalVolume() == tracker) {
207  backscattering = true;
208  } else if (postStep->GetPhysicalVolume() == cmse) {
209  // simple protection to avoid possible steps from calo towards the outer part of the detector, if allowed by geometry
210  // to be removed as soon as tracker-calo boundary becomes again the default
211  if (preStep->GetPosition().mag2() > postStep->GetPosition().mag2()) {
212  backscattering = true;
213  }
214  }
215  // store transition calo -> cmse to tag backscattering
216  if (backscattering) {
217  TrackInformation* trkinfo = static_cast<TrackInformation*>(theTrack->GetUserInformation());
218  if (!trkinfo->isInTrkFromBackscattering()) {
219  trkinfo->setInTrkFromBackscattering();
220 #ifdef EDM_ML_DEBUG
221  LogDebug("SimG4CoreApplication")
222  << "Setting flag for backscattering from CALO " << trkinfo->isInTrkFromBackscattering();
223 #endif
224  }
225  }
226  }
227  } else {
228  theTrack->SetTrackStatus(fStopAndKill);
229  isKilled = true;
230 #ifdef EDM_ML_DEBUG
231  PrintKilledTrack(theTrack, tstat);
232 #endif
233  }
234  if (nullptr != steppingVerbose) {
235  steppingVerbose->nextStep(aStep, fpSteppingManager, isKilled);
236  }
237 }
238 
239 bool Phase2SteppingAction::isLowEnergy(const G4LogicalVolume* lv, const G4Track* theTrack) const {
240  const double ekin = theTrack->GetKineticEnergy();
241  int pCode = theTrack->GetDefinition()->GetPDGEncoding();
242 
243  for (auto& vol : ekinVolumes) {
244  if (lv == vol) {
245  for (unsigned int i = 0; i < numberPart; ++i) {
246  if (pCode == ekinPDG[i]) {
247  return (ekin <= ekinMins[i]);
248  }
249  }
250  break;
251  }
252  }
253  return false;
254 }
255 
257  const G4PhysicalVolumeStore* pvs = G4PhysicalVolumeStore::GetInstance();
258  for (auto const& pvcite : *pvs) {
259  const std::string& pvname = (std::string)(DD4hep2DDDName::namePV(pvcite->GetName(), dd4hep_));
260  if (pvname == trackerName_) {
261  tracker = pvcite;
262  } else if (pvname == caloName_) {
263  calo = pvcite;
264  } else if (pvname == btlName_) {
265  btl = pvcite;
266  } else if (pvname == cmseName_) {
267  cmse = pvcite;
268  }
269  if (tracker && calo && btl && cmse)
270  break;
271  }
272  edm::LogVerbatim("SimG4CoreApplication")
273  << "Phase2SteppingAction: pointer for Tracker " << tracker << " and for Calo " << calo << " and for BTL " << btl;
274 
275  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
276  if (numberEkins > 0) {
277  ekinVolumes.resize(numberEkins, nullptr);
278  for (auto const& lvcite : *lvs) {
279  std::string lvname = (std::string)(DD4hep2DDDName::nameMatterLV(lvcite->GetName(), dd4hep_));
280  for (unsigned int i = 0; i < numberEkins; ++i) {
281  if (lvname == ekinNames[i]) {
282  ekinVolumes[i] = lvcite;
283  break;
284  }
285  }
286  }
287  for (unsigned int i = 0; i < numberEkins; ++i) {
288  edm::LogVerbatim("SimG4CoreApplication") << ekinVolumes[i]->GetName() << " with pointer " << ekinVolumes[i];
289  }
290  }
291 
292  if (numberPart > 0) {
293  G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
294  ekinPDG.resize(numberPart, 0);
295  for (unsigned int i = 0; i < numberPart; ++i) {
296  const G4ParticleDefinition* part = theParticleTable->FindParticle(ekinParticles[i]);
297  if (nullptr != part)
298  ekinPDG[i] = part->GetPDGEncoding();
299  edm::LogVerbatim("SimG4CoreApplication") << "Particle " << ekinParticles[i] << " with PDG code " << ekinPDG[i]
300  << " and KE cut off " << ekinMins[i] / MeV << " MeV";
301  }
302  }
303 
304  const G4RegionStore* rs = G4RegionStore::GetInstance();
305  if (numberTimes > 0) {
306  maxTimeRegions.resize(numberTimes, nullptr);
307  for (auto const& rcite : *rs) {
308  const G4String& rname = rcite->GetName();
309  for (unsigned int i = 0; i < numberTimes; ++i) {
310  if (rname == (G4String)(maxTimeNames[i])) {
311  maxTimeRegions[i] = rcite;
312  break;
313  }
314  }
315  }
316  }
317  if (ndeadRegions > 0) {
318  deadRegions.resize(ndeadRegions, nullptr);
319  for (auto const& rcite : *rs) {
320  const G4String& rname = rcite->GetName();
321  for (unsigned int i = 0; i < ndeadRegions; ++i) {
322  if (rname == (G4String)(deadRegionNames[i])) {
323  deadRegions[i] = rcite;
324  break;
325  }
326  }
327  }
328  }
329  return true;
330 }
331 
332 void Phase2SteppingAction::PrintKilledTrack(const G4Track* aTrack, const TrackStatus& tst) const {
333  std::string vname = "";
334  std::string rname = "";
335  std::string typ = " ";
336  switch (tst) {
337  case sDeadRegion:
338  typ = " in dead region ";
339  break;
340  case sOutOfTime:
341  typ = " out of time window ";
342  break;
343  case sLowEnergy:
344  typ = " low energy limit ";
345  break;
346  case sLowEnergyInVacuum:
347  typ = " low energy limit in vacuum ";
348  break;
349  case sEnergyDepNaN:
350  typ = " energy deposition is NaN ";
351  break;
352  case sVeryForward:
353  typ = " very forward track ";
354  break;
355  case sNumberOfSteps:
356  typ = " too many steps ";
357  break;
358  default:
359  break;
360  }
361  G4VPhysicalVolume* pv = aTrack->GetNextVolume();
362  vname = pv->GetLogicalVolume()->GetName();
363  rname = pv->GetLogicalVolume()->GetRegion()->GetName();
364 
365  const double ekin = aTrack->GetKineticEnergy();
366  if (ekin < 2 * CLHEP::MeV) {
367  return;
368  }
369  edm::LogWarning("SimG4CoreApplication")
370  << "Track #" << aTrack->GetTrackID() << " StepN= " << aTrack->GetCurrentStepNumber() << " "
371  << aTrack->GetDefinition()->GetParticleName() << " E(MeV)=" << ekin / CLHEP::MeV
372  << " T(ns)=" << aTrack->GetGlobalTime() / CLHEP::ns << " is killed due to " << typ << "\n LV: " << vname << " ("
373  << rname << ") at " << aTrack->GetPosition() << " step(cm)=" << aTrack->GetStep()->GetStepLength() / CLHEP::cm;
374 }
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)