CMS 3D CMS Logo

SteppingAction.cc
Go to the documentation of this file.
4 
5 #include "G4LogicalVolumeStore.hh"
6 #include "G4ParticleTable.hh"
7 #include "G4PhysicalVolumeStore.hh"
8 #include "G4RegionStore.hh"
9 #include "G4UnitsTable.hh"
10 #include "G4SystemOfUnits.hh"
11 
14 
15 //#define DebugLog
16 
18  : eventAction_(e),
19  tracker(nullptr),
20  calo(nullptr),
21  steppingVerbose(sv),
22  nWarnings(0),
23  initialized(false),
24  killBeamPipe(false),
25  hasWatcher(hasW) {
26  theCriticalEnergyForVacuum = (p.getParameter<double>("CriticalEnergyForVacuum") * CLHEP::MeV);
27  if (0.0 < theCriticalEnergyForVacuum) {
28  killBeamPipe = true;
29  }
30  theCriticalDensity = (p.getParameter<double>("CriticalDensity") * CLHEP::g / CLHEP::cm3);
31  maxZCentralCMS = p.getParameter<double>("MaxZCentralCMS") * CLHEP::m;
32  maxTrackTime = p.getParameter<double>("MaxTrackTime") * CLHEP::ns;
33  maxTrackTimeForward = p.getParameter<double>("MaxTrackTimeForward") * CLHEP::ns;
34  maxTrackTimes = p.getParameter<std::vector<double> >("MaxTrackTimes");
35  maxTimeNames = p.getParameter<std::vector<std::string> >("MaxTimeNames");
36  deadRegionNames = p.getParameter<std::vector<std::string> >("DeadRegions");
37  maxNumberOfSteps = p.getParameter<int>("MaxNumberOfSteps");
38  ekinMins = p.getParameter<std::vector<double> >("EkinThresholds");
39  ekinNames = p.getParameter<std::vector<std::string> >("EkinNames");
40  ekinParticles = p.getParameter<std::vector<std::string> >("EkinParticles");
41 
42  edm::LogVerbatim("SimG4CoreApplication")
43  << "SteppingAction:: 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;
50 
51  numberTimes = maxTrackTimes.size();
52  if (numberTimes > 0) {
53  for (unsigned int i = 0; i < numberTimes; i++) {
54  edm::LogVerbatim("SimG4CoreApplication")
55  << "SteppingAction::MaxTrackTime for " << maxTimeNames[i] << " is " << maxTrackTimes[i] << " ns ";
56  maxTrackTimes[i] *= ns;
57  }
58  }
59 
61  if (ndeadRegions > 0) {
62  edm::LogVerbatim("SimG4CoreApplication")
63  << "SteppingAction: Number of DeadRegions where all trackes are killed " << ndeadRegions;
64  for (unsigned int i = 0; i < ndeadRegions; ++i) {
65  edm::LogVerbatim("SimG4CoreApplication") << "SteppingAction: DeadRegion " << i << ". " << deadRegionNames[i];
66  }
67  }
68  numberEkins = ekinNames.size();
69  numberPart = ekinParticles.size();
70  if (0 == numberPart) {
71  numberEkins = 0;
72  }
73 
74  if (numberEkins > 0) {
75  edm::LogVerbatim("SimG4CoreApplication")
76  << "SteppingAction::Kill following " << numberPart << " particles in " << numberEkins << " volumes";
77  for (unsigned int i = 0; i < numberPart; ++i) {
78  edm::LogVerbatim("SimG4CoreApplication")
79  << "SteppingAction::Particle " << i << " " << ekinParticles[i] << " Threshold = " << ekinMins[i] << " MeV";
80  ekinMins[i] *= CLHEP::MeV;
81  }
82  for (unsigned int i = 0; i < numberEkins; ++i) {
83  edm::LogVerbatim("SimG4CoreApplication") << "SteppingAction::LogVolume[" << i << "] = " << ekinNames[i];
84  }
85  }
86 }
87 
89 
90 void SteppingAction::UserSteppingAction(const G4Step* aStep) {
91  if (!initialized) {
93  }
94 
95  //if(hasWatcher) { m_g4StepSignal(aStep); }
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  << "SteppingAction::UserSteppingAction: 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  // NaN energy deposit
115  if (edm::isNotFinite(aStep->GetTotalEnergyDeposit())) {
116  tstat = sEnergyDepNaN;
117  if (nWarnings < 5) {
118  ++nWarnings;
119  edm::LogWarning("SimG4CoreApplication")
120  << "Track #" << theTrack->GetTrackID() << " " << theTrack->GetDefinition()->GetParticleName()
121  << " E(MeV)= " << preStep->GetKineticEnergy() / MeV << " Nstep= " << theTrack->GetCurrentStepNumber()
122  << " is killed due to edep=NaN inside PV: " << preStep->GetPhysicalVolume()->GetName() << " at "
123  << theTrack->GetPosition() << " StepLen(mm)= " << aStep->GetStepLength();
124  }
125  }
126 
127  // the track is killed by the process
128  if (tstat == sKilledByProcess) {
129  if (nullptr != steppingVerbose) {
130  steppingVerbose->NextStep(aStep, fpSteppingManager, false);
131  }
132  return;
133  }
134 
135  if (sAlive == tstat && theTrack->GetCurrentStepNumber() > maxNumberOfSteps) {
136  tstat = sNumberOfSteps;
137  if (nWarnings < 5) {
138  ++nWarnings;
139  edm::LogWarning("SimG4CoreApplication")
140  << "Track #" << theTrack->GetTrackID() << " " << theTrack->GetDefinition()->GetParticleName()
141  << " E(MeV)= " << preStep->GetKineticEnergy() / MeV << " Nstep= " << theTrack->GetCurrentStepNumber()
142  << " is killed due to limit on number of steps;/n PV: " << preStep->GetPhysicalVolume()->GetName() << " at "
143  << theTrack->GetPosition() << " StepLen(mm)= " << aStep->GetStepLength();
144  }
145  }
146  const double time = theTrack->GetGlobalTime();
147 
148  // check Z-coordinate
149  if (sAlive == tstat && std::abs(theTrack->GetPosition().z()) >= maxZCentralCMS) {
151  }
152 
153  // check G4Region
154  if (sAlive == tstat) {
155  // next logical volume and next region
156  const G4LogicalVolume* lv = postStep->GetPhysicalVolume()->GetLogicalVolume();
157  const G4Region* theRegion = lv->GetRegion();
158 
159  // kill in dead regions
160  if (isInsideDeadRegion(theRegion))
161  tstat = sDeadRegion;
162 
163  // kill out of time
164  if (sAlive == tstat) {
165  if (isOutOfTimeWindow(theRegion, time))
166  tstat = sOutOfTime;
167  }
168 
169  // kill low-energy in volumes on demand
170  if (sAlive == tstat && numberEkins > 0) {
171  if (isLowEnergy(lv, theTrack))
172  tstat = sLowEnergy;
173  }
174 
175  // kill low-energy in vacuum
176  if (sAlive == tstat && killBeamPipe) {
177  if (theTrack->GetKineticEnergy() < theCriticalEnergyForVacuum &&
178  theTrack->GetDefinition()->GetPDGCharge() != 0.0 && lv->GetMaterial()->GetDensity() <= theCriticalDensity) {
179  tstat = sLowEnergyInVacuum;
180  }
181  }
182  }
183  // check transition tracker/calo
184  bool isKilled = false;
185  if (sAlive == tstat || sVeryForward == tstat) {
186  if (preStep->GetPhysicalVolume() == tracker && postStep->GetPhysicalVolume() == calo) {
187  math::XYZVectorD pos((postStep->GetPosition()).x(), (postStep->GetPosition()).y(), (postStep->GetPosition()).z());
188 
189  math::XYZTLorentzVectorD mom((postStep->GetMomentum()).x(),
190  (postStep->GetMomentum()).y(),
191  (postStep->GetMomentum()).z(),
192  postStep->GetTotalEnergy());
193 
194  uint32_t id = theTrack->GetTrackID();
195 
196  std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> p(pos, mom);
198  }
199  } else {
200  theTrack->SetTrackStatus(fStopAndKill);
201  isKilled = true;
202 #ifdef DebugLog
203  PrintKilledTrack(theTrack, tstat);
204 #endif
205  }
206  if (nullptr != steppingVerbose) {
207  steppingVerbose->NextStep(aStep, fpSteppingManager, isKilled);
208  }
209 }
210 
211 bool SteppingAction::isLowEnergy(const G4LogicalVolume* lv, const G4Track* theTrack) const {
212  const double ekin = theTrack->GetKineticEnergy();
213  int pCode = theTrack->GetDefinition()->GetPDGEncoding();
214 
215  for (auto& vol : ekinVolumes) {
216  if (lv == vol) {
217  for (unsigned int i = 0; i < numberPart; ++i) {
218  if (pCode == ekinPDG[i]) {
219  return (ekin <= ekinMins[i]);
220  }
221  }
222  break;
223  }
224  }
225  return false;
226 }
227 
229  const G4PhysicalVolumeStore* pvs = G4PhysicalVolumeStore::GetInstance();
230  for (auto const& pvcite : *pvs) {
231  const G4String& pvname = pvcite->GetName();
232  if (pvname == "Tracker" || pvname == "tracker:Tracker_1")
233  tracker = pvcite;
234  else if (pvname == "CALO" || pvname == "caloBase:CALO_1")
235  calo = pvcite;
236 
237  if (tracker && calo)
238  break;
239  }
240  edm::LogVerbatim("SimG4CoreApplication")
241  << "SteppingAction: pointer for Tracker " << tracker << " and for Calo " << calo;
242 
243  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
244  if (numberEkins > 0) {
245  ekinVolumes.resize(numberEkins, nullptr);
246  for (auto const& lvcite : *lvs) {
247  const G4String& lvname = lvcite->GetName();
248  for (unsigned int i = 0; i < numberEkins; ++i) {
249  if (lvname == (G4String)(ekinNames[i])) {
250  ekinVolumes[i] = lvcite;
251  break;
252  }
253  }
254  }
255  for (unsigned int i = 0; i < numberEkins; ++i) {
256  edm::LogVerbatim("SimG4CoreApplication") << ekinVolumes[i]->GetName() << " with pointer " << ekinVolumes[i];
257  }
258  }
259 
260  if (numberPart > 0) {
261  G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
262  ekinPDG.resize(numberPart, 0);
263  for (unsigned int i = 0; i < numberPart; ++i) {
264  const G4ParticleDefinition* part = theParticleTable->FindParticle(ekinParticles[i]);
265  if (nullptr != part)
266  ekinPDG[i] = part->GetPDGEncoding();
267  edm::LogVerbatim("SimG4CoreApplication") << "Particle " << ekinParticles[i] << " with PDG code " << ekinPDG[i]
268  << " and KE cut off " << ekinMins[i] / MeV << " MeV";
269  }
270  }
271 
272  const G4RegionStore* rs = G4RegionStore::GetInstance();
273  if (numberTimes > 0) {
274  maxTimeRegions.resize(numberTimes, nullptr);
275  for (auto const& rcite : *rs) {
276  const G4String& rname = rcite->GetName();
277  for (unsigned int i = 0; i < numberTimes; ++i) {
278  if (rname == (G4String)(maxTimeNames[i])) {
279  maxTimeRegions[i] = rcite;
280  break;
281  }
282  }
283  }
284  }
285  if (ndeadRegions > 0) {
286  deadRegions.resize(ndeadRegions, nullptr);
287  for (auto const& rcite : *rs) {
288  const G4String& rname = rcite->GetName();
289  for (unsigned int i = 0; i < ndeadRegions; ++i) {
290  if (rname == (G4String)(deadRegionNames[i])) {
291  deadRegions[i] = rcite;
292  break;
293  }
294  }
295  }
296  }
297  return true;
298 }
299 
300 void SteppingAction::PrintKilledTrack(const G4Track* aTrack, const TrackStatus& tst) const {
301  std::string vname = "";
302  std::string rname = "";
303  std::string typ = " ";
304  switch (tst) {
305  case sDeadRegion:
306  typ = " in dead region ";
307  break;
308  case sOutOfTime:
309  typ = " out of time window ";
310  break;
311  case sLowEnergy:
312  typ = " low energy limit ";
313  break;
314  case sLowEnergyInVacuum:
315  typ = " low energy limit in vacuum ";
316  break;
317  case sEnergyDepNaN:
318  typ = " energy deposition is NaN ";
319  break;
320  case sVeryForward:
321  typ = " very forward track ";
322  break;
323  case sNumberOfSteps:
324  typ = " too many steps ";
325  break;
326  default:
327  break;
328  }
329  G4VPhysicalVolume* pv = aTrack->GetNextVolume();
330  vname = pv->GetLogicalVolume()->GetName();
331  rname = pv->GetLogicalVolume()->GetRegion()->GetName();
332 
333  const double ekin = aTrack->GetKineticEnergy();
334  if (ekin < 2 * CLHEP::MeV) {
335  return;
336  }
337  edm::LogWarning("SimG4CoreApplication")
338  << "Track #" << aTrack->GetTrackID() << " StepN= " << aTrack->GetCurrentStepNumber() << " "
339  << aTrack->GetDefinition()->GetParticleName() << " E(MeV)=" << ekin / CLHEP::MeV
340  << " T(ns)=" << aTrack->GetGlobalTime() / CLHEP::ns << " is killed due to " << typ << "\n LV: " << vname << " ("
341  << rname << ") at " << aTrack->GetPosition() << " step(cm)=" << aTrack->GetStep()->GetStepLength() / CLHEP::cm;
342 }
Log< level::Info, true > LogVerbatim
void UserSteppingAction(const G4Step *aStep) final
const G4VPhysicalVolume * tracker
std::vector< int > ekinPDG
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
std::vector< const G4Region * > deadRegions
double theCriticalDensity
SimActivityRegistry::G4StepSignal m_g4StepSignal
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
~SteppingAction() override
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
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
TrackStatus
G4int maxNumberOfSteps
void PrintKilledTrack(const G4Track *, const TrackStatus &) const
std::vector< std::string > deadRegionNames
const G4VPhysicalVolume * calo
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
double theCriticalEnergyForVacuum
EventAction * eventAction_
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 addTkCaloStateInfo(uint32_t t, const std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > &p)
Definition: EventAction.cc:71
void NextStep(const G4Step *, const G4SteppingManager *ptr, bool isKilled) const
SteppingAction(EventAction *ea, const edm::ParameterSet &ps, const CMSSteppingVerbose *, bool hasW)
Log< level::Warning, false > LogWarning
unsigned int ndeadRegions
std::vector< std::string > ekinParticles
Definition: Common.h:9
std::vector< const G4Region * > maxTimeRegions
unsigned int nWarnings