CMS 3D CMS Logo

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