CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
61  ndeadRegions = deadRegionNames.size();
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) {
141  tstat = (time > maxTrackTimeForward) ? sOutOfTime : sVeryForward;
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  if (sAlive == tstat || sVeryForward == tstat) {
176  if (isThisVolume(preStep->GetTouchable(), tracker) && isThisVolume(postStep->GetTouchable(), calo)) {
177  math::XYZVectorD pos((preStep->GetPosition()).x(), (preStep->GetPosition()).y(), (preStep->GetPosition()).z());
178 
179  math::XYZTLorentzVectorD mom((preStep->GetMomentum()).x(),
180  (preStep->GetMomentum()).y(),
181  (preStep->GetMomentum()).z(),
182  preStep->GetTotalEnergy());
183 
184  uint32_t id = theTrack->GetTrackID();
185 
186  std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> p(pos, mom);
188  }
189  } else {
190  theTrack->SetTrackStatus(fStopAndKill);
191 #ifdef DebugLog
192  PrintKilledTrack(theTrack, tstat);
193 #endif
194  }
195  if (nullptr != steppingVerbose) {
196  steppingVerbose->NextStep(aStep, fpSteppingManager, (1 < tstat));
197  }
198 }
199 
200 bool SteppingAction::isLowEnergy(const G4LogicalVolume* lv, const G4Track* theTrack) const {
201  const double ekin = theTrack->GetKineticEnergy();
202  int pCode = theTrack->GetDefinition()->GetPDGEncoding();
203 
204  for (auto& vol : ekinVolumes) {
205  if (lv == vol) {
206  for (unsigned int i = 0; i < numberPart; ++i) {
207  if (pCode == ekinPDG[i]) {
208  return (ekin <= ekinMins[i]);
209  }
210  }
211  break;
212  }
213  }
214  return false;
215 }
216 
218  const G4PhysicalVolumeStore* pvs = G4PhysicalVolumeStore::GetInstance();
219  for (auto const& pvcite : *pvs) {
220  const G4String& pvname = pvcite->GetName();
221  if (pvname == "Tracker")
222  tracker = pvcite;
223  else if (pvname == "CALO")
224  calo = pvcite;
225 
226  if (tracker && calo)
227  break;
228  }
229  edm::LogVerbatim("SimG4CoreApplication")
230  << "SteppingAction: pointer for Tracker " << tracker << " and for Calo " << calo;
231 
232  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
233  if (numberEkins > 0) {
234  ekinVolumes.resize(numberEkins, nullptr);
235  for (auto const& lvcite : *lvs) {
236  const G4String& lvname = lvcite->GetName();
237  for (unsigned int i = 0; i < numberEkins; ++i) {
238  if (lvname == (G4String)(ekinNames[i])) {
239  ekinVolumes[i] = lvcite;
240  break;
241  }
242  }
243  }
244  for (unsigned int i = 0; i < numberEkins; ++i) {
245  edm::LogVerbatim("SimG4CoreApplication") << ekinVolumes[i]->GetName() << " with pointer " << ekinVolumes[i];
246  }
247  }
248 
249  if (numberPart > 0) {
250  G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
251  ekinPDG.resize(numberPart, 0);
252  for (unsigned int i = 0; i < numberPart; ++i) {
253  const G4ParticleDefinition* part = theParticleTable->FindParticle(ekinParticles[i]);
254  if (nullptr != part)
255  ekinPDG[i] = part->GetPDGEncoding();
256  edm::LogVerbatim("SimG4CoreApplication") << "Particle " << ekinParticles[i] << " with PDG code " << ekinPDG[i]
257  << " and KE cut off " << ekinMins[i] / MeV << " MeV";
258  }
259  }
260 
261  const G4RegionStore* rs = G4RegionStore::GetInstance();
262  if (numberTimes > 0) {
263  maxTimeRegions.resize(numberTimes, nullptr);
264  for (auto const& rcite : *rs) {
265  const G4String& rname = rcite->GetName();
266  for (unsigned int i = 0; i < numberTimes; ++i) {
267  if (rname == (G4String)(maxTimeNames[i])) {
268  maxTimeRegions[i] = rcite;
269  break;
270  }
271  }
272  }
273  }
274  if (ndeadRegions > 0) {
275  deadRegions.resize(ndeadRegions, nullptr);
276  for (auto const& rcite : *rs) {
277  const G4String& rname = rcite->GetName();
278  for (unsigned int i = 0; i < ndeadRegions; ++i) {
279  if (rname == (G4String)(deadRegionNames[i])) {
280  deadRegions[i] = rcite;
281  break;
282  }
283  }
284  }
285  }
286  return true;
287 }
288 
289 void SteppingAction::PrintKilledTrack(const G4Track* aTrack, const TrackStatus& tst) const {
290  std::string vname = "";
291  std::string rname = "";
292  std::string typ = " ";
293  switch (tst) {
294  case sDeadRegion:
295  typ = " in dead region ";
296  break;
297  case sOutOfTime:
298  typ = " out of time window ";
299  break;
300  case sLowEnergy:
301  typ = " low energy limit ";
302  break;
303  case sLowEnergyInVacuum:
304  typ = " low energy limit in vacuum ";
305  break;
306  case sEnergyDepNaN:
307  typ = " energy deposition is NaN ";
308  break;
309  case sVeryForward:
310  typ = " very forward track ";
311  break;
312  case sNumberOfSteps:
313  typ = " too many steps ";
314  break;
315  default:
316  break;
317  }
318  G4VPhysicalVolume* pv = aTrack->GetNextVolume();
319  vname = pv->GetLogicalVolume()->GetName();
320  rname = pv->GetLogicalVolume()->GetRegion()->GetName();
321 
322  const double ekin = aTrack->GetKineticEnergy();
323  if (ekin < 2 * CLHEP::MeV) {
324  return;
325  }
326  edm::LogWarning("SimG4CoreApplication")
327  << "Track #" << aTrack->GetTrackID() << " StepN= " << aTrack->GetCurrentStepNumber() << " "
328  << aTrack->GetDefinition()->GetParticleName() << " E(MeV)=" << ekin / CLHEP::MeV
329  << " T(ns)=" << aTrack->GetGlobalTime() / CLHEP::ns << " is killed due to " << typ << "\n LV: " << vname << " ("
330  << rname << ") at " << aTrack->GetPosition() << " step(cm)=" << aTrack->GetStep()->GetStepLength() / CLHEP::cm;
331 }
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
void NextStep(const G4Step *, const G4SteppingManager *ptr, bool isKilled) const
unsigned int numberPart
const CMSSteppingVerbose * steppingVerbose
unsigned int numberEkins
bool isLowEnergy(const G4LogicalVolume *, const G4Track *) 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
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
TrackStatus
bool isOutOfTimeWindow(const G4Region *reg, const double &time) const
G4int maxNumberOfSteps
std::vector< std::string > deadRegionNames
const G4VPhysicalVolume * calo
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int numberTimes
std::vector< double > maxTrackTimes
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
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double maxTrackTimeForward
const G4String rname[NREG]
void addTkCaloStateInfo(uint32_t t, const std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > &p)
Definition: EventAction.cc:72
bool isInsideDeadRegion(const G4Region *reg) const
SteppingAction(EventAction *ea, const edm::ParameterSet &ps, const CMSSteppingVerbose *, bool hasW)
Log< level::Warning, false > LogWarning
void PrintKilledTrack(const G4Track *, const TrackStatus &) const
unsigned int ndeadRegions
std::vector< std::string > ekinParticles
bool isThisVolume(const G4VTouchable *touch, const G4VPhysicalVolume *pv) const
std::vector< const G4Region * > maxTimeRegions
unsigned int nWarnings