CMS 3D CMS Logo

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