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