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  trkinfo->setIdAtBTLentrance(theTrack->GetTrackID());
172 #ifdef DebugLog
173  LogDebug("SimG4CoreApplication") << "Setting flag for Tracker -> BTL " << trkinfo->isFromTtoBTL()
174  << " IdAtBTLentrance = " << trkinfo->idAtBTLentrance();
175 #endif
176  } else {
177  trkinfo->setBTLlooper();
178  trkinfo->setIdAtBTLentrance(theTrack->GetTrackID());
179 #ifdef DebugLog
180  LogDebug("SimG4CoreApplication") << "Setting flag for BTL looper " << trkinfo->isBTLlooper();
181  trkinfo->Print();
182 #endif
183  }
184  } else if (preStep->GetPhysicalVolume() == btl && postStep->GetPhysicalVolume() == tracker) {
185  // store transition BTL -> tracker
186  TrackInformation* trkinfo = static_cast<TrackInformation*>(theTrack->GetUserInformation());
187  if (!trkinfo->isFromBTLtoT()) {
188  trkinfo->setFromBTLtoT();
189 #ifdef DebugLog
190  LogDebug("SimG4CoreApplication") << "Setting flag for BTL -> Tracker " << trkinfo->isFromBTLtoT();
191 #endif
192  }
193  } else if (preStep->GetPhysicalVolume() == tracker && postStep->GetPhysicalVolume() == calo) {
194  // store transition tracker -> calo
195  TrackInformation* trkinfo = static_cast<TrackInformation*>(theTrack->GetUserInformation());
196  if (!trkinfo->crossedBoundary()) {
197  trkinfo->setCrossedBoundary(theTrack);
198  }
199  }
200  } else {
201  theTrack->SetTrackStatus(fStopAndKill);
202  isKilled = true;
203 #ifdef DebugLog
204  PrintKilledTrack(theTrack, tstat);
205 #endif
206  }
207  if (nullptr != steppingVerbose) {
208  steppingVerbose->nextStep(aStep, fpSteppingManager, isKilled);
209  }
210 }
211 
212 bool Phase2SteppingAction::isLowEnergy(const G4LogicalVolume* lv, const G4Track* theTrack) const {
213  const double ekin = theTrack->GetKineticEnergy();
214  int pCode = theTrack->GetDefinition()->GetPDGEncoding();
215 
216  for (auto& vol : ekinVolumes) {
217  if (lv == vol) {
218  for (unsigned int i = 0; i < numberPart; ++i) {
219  if (pCode == ekinPDG[i]) {
220  return (ekin <= ekinMins[i]);
221  }
222  }
223  break;
224  }
225  }
226  return false;
227 }
228 
230  const G4PhysicalVolumeStore* pvs = G4PhysicalVolumeStore::GetInstance();
231  for (auto const& pvcite : *pvs) {
232  const G4String& pvname = pvcite->GetName();
233  if (pvname == "Tracker" || pvname == "tracker:Tracker_1") {
234  tracker = pvcite;
235  } else if (pvname == "CALO" || pvname == "caloBase:CALO_1") {
236  calo = pvcite;
237  } else if (pvname == "BarrelTimingLayer" || pvname == "btl:BarrelTimingLayer_1") {
238  btl = pvcite;
239  }
240  if (tracker && calo && btl)
241  break;
242  }
243  edm::LogVerbatim("SimG4CoreApplication")
244  << "Phase2SteppingAction: pointer for Tracker " << tracker << " and for Calo " << calo << " and for BTL " << btl;
245 
246  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
247  if (numberEkins > 0) {
248  ekinVolumes.resize(numberEkins, nullptr);
249  for (auto const& lvcite : *lvs) {
250  const G4String& lvname = lvcite->GetName();
251  for (unsigned int i = 0; i < numberEkins; ++i) {
252  if (lvname == (G4String)(ekinNames[i])) {
253  ekinVolumes[i] = lvcite;
254  break;
255  }
256  }
257  }
258  for (unsigned int i = 0; i < numberEkins; ++i) {
259  edm::LogVerbatim("SimG4CoreApplication") << ekinVolumes[i]->GetName() << " with pointer " << ekinVolumes[i];
260  }
261  }
262 
263  if (numberPart > 0) {
264  G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
265  ekinPDG.resize(numberPart, 0);
266  for (unsigned int i = 0; i < numberPart; ++i) {
267  const G4ParticleDefinition* part = theParticleTable->FindParticle(ekinParticles[i]);
268  if (nullptr != part)
269  ekinPDG[i] = part->GetPDGEncoding();
270  edm::LogVerbatim("SimG4CoreApplication") << "Particle " << ekinParticles[i] << " with PDG code " << ekinPDG[i]
271  << " and KE cut off " << ekinMins[i] / MeV << " MeV";
272  }
273  }
274 
275  const G4RegionStore* rs = G4RegionStore::GetInstance();
276  if (numberTimes > 0) {
277  maxTimeRegions.resize(numberTimes, nullptr);
278  for (auto const& rcite : *rs) {
279  const G4String& rname = rcite->GetName();
280  for (unsigned int i = 0; i < numberTimes; ++i) {
281  if (rname == (G4String)(maxTimeNames[i])) {
282  maxTimeRegions[i] = rcite;
283  break;
284  }
285  }
286  }
287  }
288  if (ndeadRegions > 0) {
289  deadRegions.resize(ndeadRegions, nullptr);
290  for (auto const& rcite : *rs) {
291  const G4String& rname = rcite->GetName();
292  for (unsigned int i = 0; i < ndeadRegions; ++i) {
293  if (rname == (G4String)(deadRegionNames[i])) {
294  deadRegions[i] = rcite;
295  break;
296  }
297  }
298  }
299  }
300  return true;
301 }
302 
303 void Phase2SteppingAction::PrintKilledTrack(const G4Track* aTrack, const TrackStatus& tst) const {
304  std::string vname = "";
305  std::string rname = "";
306  std::string typ = " ";
307  switch (tst) {
308  case sDeadRegion:
309  typ = " in dead region ";
310  break;
311  case sOutOfTime:
312  typ = " out of time window ";
313  break;
314  case sLowEnergy:
315  typ = " low energy limit ";
316  break;
317  case sLowEnergyInVacuum:
318  typ = " low energy limit in vacuum ";
319  break;
320  case sEnergyDepNaN:
321  typ = " energy deposition is NaN ";
322  break;
323  case sVeryForward:
324  typ = " very forward track ";
325  break;
326  case sNumberOfSteps:
327  typ = " too many steps ";
328  break;
329  default:
330  break;
331  }
332  G4VPhysicalVolume* pv = aTrack->GetNextVolume();
333  vname = pv->GetLogicalVolume()->GetName();
334  rname = pv->GetLogicalVolume()->GetRegion()->GetName();
335 
336  const double ekin = aTrack->GetKineticEnergy();
337  if (ekin < 2 * CLHEP::MeV) {
338  return;
339  }
340  edm::LogWarning("SimG4CoreApplication")
341  << "Track #" << aTrack->GetTrackID() << " StepN= " << aTrack->GetCurrentStepNumber() << " "
342  << aTrack->GetDefinition()->GetParticleName() << " E(MeV)=" << ekin / CLHEP::MeV
343  << " T(ns)=" << aTrack->GetGlobalTime() / CLHEP::ns << " is killed due to " << typ << "\n LV: " << vname << " ("
344  << rname << ") at " << aTrack->GetPosition() << " step(cm)=" << aTrack->GetStep()->GetStepLength() / CLHEP::cm;
345 }
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
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
void setIdAtBTLentrance(int id)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< G4LogicalVolume * > ekinVolumes
const G4VPhysicalVolume * btl
bool isFromBTLtoT() const
void Print() const override
std::vector< std::string > ekinParticles
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
int idAtBTLentrance() 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
std::vector< double > maxTrackTimes
void UserSteppingAction(const G4Step *aStep) final
#define LogDebug(id)