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  : trackManager_(stm), 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  math::XYZVectorD pos((postStep->GetPosition()).x(), (postStep->GetPosition()).y(), (postStep->GetPosition()).z());
166 
167  math::XYZTLorentzVectorD mom((postStep->GetMomentum()).x(),
168  (postStep->GetMomentum()).y(),
169  (postStep->GetMomentum()).z(),
170  postStep->GetTotalEnergy());
171 
172  // record intersection
173  uint32_t id = theTrack->GetTrackID();
174  std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> p(pos, mom);
176  }
177  } else {
178  theTrack->SetTrackStatus(fStopAndKill);
179  isKilled = true;
180 #ifdef DebugLog
181  PrintKilledTrack(theTrack, tstat);
182 #endif
183  }
184  if (nullptr != steppingVerbose) {
185  steppingVerbose->nextStep(aStep, fpSteppingManager, isKilled);
186  }
187 }
188 
189 bool SteppingAction::isLowEnergy(const G4LogicalVolume* lv, const G4Track* theTrack) const {
190  const double ekin = theTrack->GetKineticEnergy();
191  int pCode = theTrack->GetDefinition()->GetPDGEncoding();
192 
193  for (auto& vol : ekinVolumes) {
194  if (lv == vol) {
195  for (unsigned int i = 0; i < numberPart; ++i) {
196  if (pCode == ekinPDG[i]) {
197  return (ekin <= ekinMins[i]);
198  }
199  }
200  break;
201  }
202  }
203  return false;
204 }
205 
207  const G4PhysicalVolumeStore* pvs = G4PhysicalVolumeStore::GetInstance();
208  for (auto const& pvcite : *pvs) {
209  const G4String& pvname = pvcite->GetName();
210  if (pvname == "Tracker" || pvname == "tracker:Tracker_1") {
211  tracker = pvcite;
212  } else if (pvname == "CALO" || pvname == "caloBase:CALO_1") {
213  calo = pvcite;
214  }
215  if (tracker && calo)
216  break;
217  }
218  edm::LogVerbatim("SimG4CoreApplication")
219  << "SteppingAction: pointer for Tracker " << tracker << " and for Calo " << calo;
220 
221  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
222  if (numberEkins > 0) {
223  ekinVolumes.resize(numberEkins, nullptr);
224  for (auto const& lvcite : *lvs) {
225  const G4String& lvname = lvcite->GetName();
226  for (unsigned int i = 0; i < numberEkins; ++i) {
227  if (lvname == (G4String)(ekinNames[i])) {
228  ekinVolumes[i] = lvcite;
229  break;
230  }
231  }
232  }
233  for (unsigned int i = 0; i < numberEkins; ++i) {
234  edm::LogVerbatim("SimG4CoreApplication") << ekinVolumes[i]->GetName() << " with pointer " << ekinVolumes[i];
235  }
236  }
237 
238  if (numberPart > 0) {
239  G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
240  ekinPDG.resize(numberPart, 0);
241  for (unsigned int i = 0; i < numberPart; ++i) {
242  const G4ParticleDefinition* part = theParticleTable->FindParticle(ekinParticles[i]);
243  if (nullptr != part)
244  ekinPDG[i] = part->GetPDGEncoding();
245  edm::LogVerbatim("SimG4CoreApplication") << "Particle " << ekinParticles[i] << " with PDG code " << ekinPDG[i]
246  << " and KE cut off " << ekinMins[i] / MeV << " MeV";
247  }
248  }
249 
250  const G4RegionStore* rs = G4RegionStore::GetInstance();
251  if (numberTimes > 0) {
252  maxTimeRegions.resize(numberTimes, nullptr);
253  for (auto const& rcite : *rs) {
254  const G4String& rname = rcite->GetName();
255  for (unsigned int i = 0; i < numberTimes; ++i) {
256  if (rname == (G4String)(maxTimeNames[i])) {
257  maxTimeRegions[i] = rcite;
258  break;
259  }
260  }
261  }
262  }
263  if (ndeadRegions > 0) {
264  deadRegions.resize(ndeadRegions, nullptr);
265  for (auto const& rcite : *rs) {
266  const G4String& rname = rcite->GetName();
267  for (unsigned int i = 0; i < ndeadRegions; ++i) {
268  if (rname == (G4String)(deadRegionNames[i])) {
269  deadRegions[i] = rcite;
270  break;
271  }
272  }
273  }
274  }
275  return true;
276 }
277 
278 void SteppingAction::PrintKilledTrack(const G4Track* aTrack, const TrackStatus& tst) const {
279  std::string vname = "";
280  std::string rname = "";
281  std::string typ = " ";
282  switch (tst) {
283  case sDeadRegion:
284  typ = " in dead region ";
285  break;
286  case sOutOfTime:
287  typ = " out of time window ";
288  break;
289  case sLowEnergy:
290  typ = " low energy limit ";
291  break;
292  case sLowEnergyInVacuum:
293  typ = " low energy limit in vacuum ";
294  break;
295  case sEnergyDepNaN:
296  typ = " energy deposition is NaN ";
297  break;
298  case sVeryForward:
299  typ = " very forward track ";
300  break;
301  case sNumberOfSteps:
302  typ = " too many steps ";
303  break;
304  default:
305  break;
306  }
307  G4VPhysicalVolume* pv = aTrack->GetNextVolume();
308  vname = pv->GetLogicalVolume()->GetName();
309  rname = pv->GetLogicalVolume()->GetRegion()->GetName();
310 
311  const double ekin = aTrack->GetKineticEnergy();
312  if (ekin < 2 * CLHEP::MeV) {
313  return;
314  }
315  edm::LogWarning("SimG4CoreApplication")
316  << "Track #" << aTrack->GetTrackID() << " StepN= " << aTrack->GetCurrentStepNumber() << " "
317  << aTrack->GetDefinition()->GetParticleName() << " E(MeV)=" << ekin / CLHEP::MeV
318  << " T(ns)=" << aTrack->GetGlobalTime() / CLHEP::ns << " is killed due to " << typ << "\n LV: " << vname << " ("
319  << rname << ") at " << aTrack->GetPosition() << " step(cm)=" << aTrack->GetStep()->GetStepLength() / CLHEP::cm;
320 }
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
SteppingAction(SimTrackManager *, const CMSSteppingVerbose *, const edm::ParameterSet &, bool hasW)
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
void addTkCaloStateInfo(uint32_t t, const std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > &p)
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
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
SimTrackManager * trackManager_
Definition: Common.h:9
std::vector< const G4Region * > maxTimeRegions
unsigned int nWarnings