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),
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  maxTrackTime = p.getParameter<double>("MaxTrackTime") * CLHEP::ns;
33  maxTrackTimes = p.getParameter<std::vector<double> >("MaxTrackTimes");
34  maxTimeNames = p.getParameter<std::vector<std::string> >("MaxTimeNames");
35  deadRegionNames = p.getParameter<std::vector<std::string> >("DeadRegions");
36  ekinMins = p.getParameter<std::vector<double> >("EkinThresholds");
37  ekinNames = p.getParameter<std::vector<std::string> >("EkinNames");
38  ekinParticles = p.getParameter<std::vector<std::string> >("EkinParticles");
39 
40  edm::LogVerbatim("SimG4CoreApplication")
41  << "SteppingAction:: KillBeamPipe = " << killBeamPipe
42  << " CriticalDensity = " << theCriticalDensity * CLHEP::cm3 / CLHEP::g << " g/cm3;"
43  << " CriticalEnergyForVacuum = " << theCriticalEnergyForVacuum / CLHEP::MeV << " Mev;"
44  << " MaxTrackTime = " << maxTrackTime / CLHEP::ns << " ns";
45 
46  numberTimes = maxTrackTimes.size();
47  if (numberTimes > 0) {
48  for (unsigned int i = 0; i < numberTimes; i++) {
49  edm::LogVerbatim("SimG4CoreApplication")
50  << "SteppingAction::MaxTrackTime for " << maxTimeNames[i] << " is " << maxTrackTimes[i] << " ns ";
51  maxTrackTimes[i] *= ns;
52  }
53  }
54 
55  ndeadRegions = deadRegionNames.size();
56  if (ndeadRegions > 0) {
57  edm::LogVerbatim("SimG4CoreApplication")
58  << "SteppingAction: Number of DeadRegions where all trackes are killed " << ndeadRegions;
59  for (unsigned int i = 0; i < ndeadRegions; ++i) {
60  edm::LogVerbatim("SimG4CoreApplication") << "SteppingAction: DeadRegion " << i << ". " << deadRegionNames[i];
61  }
62  }
63  numberEkins = ekinNames.size();
64  numberPart = ekinParticles.size();
65  if (0 == numberPart) {
66  numberEkins = 0;
67  }
68 
69  if (numberEkins > 0) {
70  edm::LogVerbatim("SimG4CoreApplication")
71  << "SteppingAction::Kill following " << numberPart << " particles in " << numberEkins << " volumes";
72  for (unsigned int i = 0; i < numberPart; ++i) {
73  edm::LogVerbatim("SimG4CoreApplication")
74  << "SteppingAction::Particle " << i << " " << ekinParticles[i] << " Threshold = " << ekinMins[i] << " MeV";
75  ekinMins[i] *= CLHEP::MeV;
76  }
77  for (unsigned int i = 0; i < numberEkins; ++i) {
78  edm::LogVerbatim("SimG4CoreApplication") << "SteppingAction::LogVolume[" << i << "] = " << ekinNames[i];
79  }
80  }
81 }
82 
84 
85 void SteppingAction::UserSteppingAction(const G4Step* aStep) {
86  if (!initialized) {
88  }
89 
90  //if(hasWatcher) { m_g4StepSignal(aStep); }
91  m_g4StepSignal(aStep);
92 
93  G4Track* theTrack = aStep->GetTrack();
94  TrackStatus tstat = (theTrack->GetTrackStatus() == fAlive) ? sAlive : sKilledByProcess;
95 
96  G4StepPoint* postStep = aStep->GetPostStepPoint();
97 
98  if (sAlive == tstat && postStep->GetPhysicalVolume() != nullptr) {
99  G4StepPoint* preStep = aStep->GetPreStepPoint();
100  const G4Region* theRegion = preStep->GetPhysicalVolume()->GetLogicalVolume()->GetRegion();
101 
102  // kill in dead regions
103  if (isInsideDeadRegion(theRegion)) {
104  tstat = sDeadRegion;
105  }
106 
107  // NaN energy deposit
108  if (sAlive == tstat && edm::isNotFinite(aStep->GetTotalEnergyDeposit())) {
109  tstat = sEnergyDepNaN;
110  if (nWarnings < 20) {
111  ++nWarnings;
112  edm::LogWarning("SimG4CoreApplication")
113  << "Track #" << theTrack->GetTrackID() << " " << theTrack->GetDefinition()->GetParticleName()
114  << " E(MeV)= " << preStep->GetKineticEnergy() / MeV
115  << " is killed due to edep=NaN inside PV: " << preStep->GetPhysicalVolume()->GetName() << " at "
116  << theTrack->GetPosition() << " StepLen(mm)= " << aStep->GetStepLength();
117  }
118  }
119 
120  // kill out of time
121  if (sAlive == tstat && isOutOfTimeWindow(theTrack, theRegion)) {
122  tstat = sOutOfTime;
123  }
124 
125  // kill low-energy in volumes on demand
126  if (sAlive == tstat && numberEkins > 0 && isLowEnergy(aStep)) {
127  tstat = sLowEnergy;
128  }
129 
130  // kill low-energy in vacuum
131  G4double kinEnergy = theTrack->GetKineticEnergy();
132  if (sAlive == tstat && killBeamPipe && kinEnergy < theCriticalEnergyForVacuum &&
133  theTrack->GetDefinition()->GetPDGCharge() != 0.0 && kinEnergy > 0.0 &&
134  theTrack->GetNextVolume()->GetLogicalVolume()->GetMaterial()->GetDensity() <= theCriticalDensity) {
135  tstat = sLowEnergyInVacuum;
136  }
137 
138  // check transition tracker/calo
139  if (sAlive == tstat) {
140  if (isThisVolume(preStep->GetTouchable(), tracker) && isThisVolume(postStep->GetTouchable(), calo)) {
141  math::XYZVectorD pos((preStep->GetPosition()).x(), (preStep->GetPosition()).y(), (preStep->GetPosition()).z());
142 
143  math::XYZTLorentzVectorD mom((preStep->GetMomentum()).x(),
144  (preStep->GetMomentum()).y(),
145  (preStep->GetMomentum()).z(),
146  preStep->GetTotalEnergy());
147 
148  uint32_t id = theTrack->GetTrackID();
149 
150  std::pair<math::XYZVectorD, math::XYZTLorentzVectorD> p(pos, mom);
152  }
153  } else {
154  theTrack->SetTrackStatus(fStopAndKill);
155 #ifdef DebugLog
156  PrintKilledTrack(theTrack, tstat);
157 #endif
158  }
159  }
160  if (nullptr != steppingVerbose) {
161  steppingVerbose->NextStep(aStep, fpSteppingManager, (1 < tstat));
162  }
163 }
164 
165 bool SteppingAction::isLowEnergy(const G4Step* aStep) const {
166  bool flag = false;
167  const G4StepPoint* sp = aStep->GetPostStepPoint();
168  G4LogicalVolume* lv = sp->GetPhysicalVolume()->GetLogicalVolume();
169  for (unsigned int i = 0; i < numberEkins; ++i) {
170  if (lv == ekinVolumes[i]) {
171  flag = true;
172  break;
173  }
174  }
175  if (flag) {
176  double ekin = sp->GetKineticEnergy();
177  int pCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
178  for (unsigned int i = 0; i < numberPart; ++i) {
179  if (pCode == ekinPDG[i]) {
180  return (ekin <= ekinMins[i]) ? true : false;
181  }
182  }
183  }
184  return false;
185 }
186 
188  const G4PhysicalVolumeStore* pvs = G4PhysicalVolumeStore::GetInstance();
189  if (pvs) {
190  std::vector<G4VPhysicalVolume*>::const_iterator pvcite;
191  for (pvcite = pvs->begin(); pvcite != pvs->end(); ++pvcite) {
192  if ((*pvcite)->GetName() == "Tracker")
193  tracker = (*pvcite);
194  if ((*pvcite)->GetName() == "CALO")
195  calo = (*pvcite);
196  if (tracker && calo)
197  break;
198  }
199  if (tracker || calo) {
200  edm::LogVerbatim("SimG4CoreApplication") << "Pointer for Tracker " << tracker << " and for Calo " << calo;
201  if (tracker)
202  LogDebug("SimG4CoreApplication") << "Tracker vol name " << tracker->GetName();
203  if (calo)
204  LogDebug("SimG4CoreApplication") << "Calorimeter vol name " << calo->GetName();
205  }
206  }
207 
208  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
209  if (numberEkins > 0) {
210  if (lvs) {
211  ekinVolumes.resize(numberEkins, nullptr);
212  std::vector<G4LogicalVolume*>::const_iterator lvcite;
213  for (lvcite = lvs->begin(); lvcite != lvs->end(); ++lvcite) {
214  for (unsigned int i = 0; i < numberEkins; ++i) {
215  if ((*lvcite)->GetName() == (G4String)(ekinNames[i])) {
216  ekinVolumes[i] = (*lvcite);
217  break;
218  }
219  }
220  }
221  }
222  for (unsigned int i = 0; i < numberEkins; ++i) {
223  edm::LogVerbatim("SimG4CoreApplication") << ekinVolumes[i]->GetName() << " with pointer " << ekinVolumes[i];
224  }
225  }
226 
227  if (numberPart > 0) {
228  G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
229  G4String partName;
230  ekinPDG.resize(numberPart, 0);
231  for (unsigned int i = 0; i < numberPart; ++i) {
232  ekinPDG[i] = theParticleTable->FindParticle(partName = ekinParticles[i])->GetPDGEncoding();
233  edm::LogVerbatim("SimG4CoreApplication") << "Particle " << ekinParticles[i] << " with PDG code " << ekinPDG[i]
234  << " and KE cut off " << ekinMins[i] / MeV << " MeV";
235  }
236  }
237 
238  const G4RegionStore* rs = G4RegionStore::GetInstance();
239  if (numberTimes > 0) {
240  maxTimeRegions.resize(numberTimes, nullptr);
241  std::vector<G4Region*>::const_iterator rcite;
242  for (rcite = rs->begin(); rcite != rs->end(); ++rcite) {
243  for (unsigned int i = 0; i < numberTimes; ++i) {
244  if ((*rcite)->GetName() == (G4String)(maxTimeNames[i])) {
245  maxTimeRegions[i] = (*rcite);
246  break;
247  }
248  }
249  }
250  }
251  if (ndeadRegions > 0) {
252  deadRegions.resize(ndeadRegions, nullptr);
253  std::vector<G4Region*>::const_iterator rcite;
254  for (rcite = rs->begin(); rcite != rs->end(); ++rcite) {
255  for (unsigned int i = 0; i < ndeadRegions; ++i) {
256  if ((*rcite)->GetName() == (G4String)(deadRegionNames[i])) {
257  deadRegions[i] = (*rcite);
258  break;
259  }
260  }
261  }
262  }
263  return true;
264 }
265 
266 void SteppingAction::PrintKilledTrack(const G4Track* aTrack, const TrackStatus& tst) const {
267  std::string vname = "";
268  std::string rname = "";
269  std::string typ = " ";
270  switch (tst) {
271  case sKilledByProcess:
272  typ = " G4Process ";
273  break;
274  case sDeadRegion:
275  typ = " in dead region ";
276  break;
277  case sOutOfTime:
278  typ = " out of time window ";
279  break;
280  case sLowEnergy:
281  typ = " low energy limit ";
282  break;
283  case sLowEnergyInVacuum:
284  typ = " low energy limit in vacuum ";
285  break;
286  case sEnergyDepNaN:
287  typ = " energy deposition is NaN ";
288  break;
289  default:
290  break;
291  }
292  G4VPhysicalVolume* pv = aTrack->GetNextVolume();
293  if (pv) {
294  vname = pv->GetLogicalVolume()->GetName();
295  rname = pv->GetLogicalVolume()->GetRegion()->GetName();
296  }
297 
298  edm::LogVerbatim("SimG4CoreApplication")
299  << "Track #" << aTrack->GetTrackID() << " " << aTrack->GetDefinition()->GetParticleName()
300  << " E(MeV)= " << aTrack->GetKineticEnergy() / MeV << " T(ns)= " << aTrack->GetGlobalTime() / ns
301  << " is killed due to " << typ << " inside LV: " << vname << " (" << rname << ") at " << aTrack->GetPosition();
302 }
#define LogDebug(id)
T getParameter(std::string const &) const
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
#define nullptr
unsigned int numberPart
const CMSSteppingVerbose * steppingVerbose
bool isOutOfTimeWindow(G4Track *theTrack, const G4Region *reg) const
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
const double MeV
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
TrackStatus
std::vector< std::string > deadRegionNames
const G4VPhysicalVolume * calo
def pv(vc)
Definition: MetAnalyzer.py:7
unsigned int numberTimes
std::vector< double > maxTrackTimes
std::vector< G4LogicalVolume * > ekinVolumes
double theCriticalEnergyForVacuum
EventAction * eventAction_
std::vector< double > ekinMins
std::vector< std::string > maxTimeNames
std::vector< std::string > ekinNames
const G4String rname[NREG]
void addTkCaloStateInfo(uint32_t t, const std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > &p)
Definition: EventAction.cc:72
bool isLowEnergy(const G4Step *aStep) const
bool isInsideDeadRegion(const G4Region *reg) const
SteppingAction(EventAction *ea, const edm::ParameterSet &ps, const CMSSteppingVerbose *, bool hasW)
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