CMS 3D CMS Logo

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