CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Public Attributes | Private Member Functions | Private Attributes
SteppingAction Class Reference

#include <SteppingAction.h>

Inheritance diagram for SteppingAction:

Public Member Functions

 SteppingAction (EventAction *ea, const edm::ParameterSet &ps)
 
virtual void UserSteppingAction (const G4Step *aStep)
 
virtual ~SteppingAction ()
 

Public Attributes

SimActivityRegistry::G4StepSignal m_g4StepSignal
 

Private Member Functions

bool catchLongLived (G4Track *theTrack, const G4Region *reg) const
 
bool initPointer ()
 
bool isThisVolume (const G4VTouchable *touch, G4VPhysicalVolume *pv) const
 
bool killInsideDeadRegion (G4Track *theTrack, const G4Region *reg) const
 
bool killLowEnergy (const G4Step *aStep) const
 
void PrintKilledTrack (const G4Track *, const std::string &) const
 

Private Attributes

G4VPhysicalVolume * calo
 
std::vector< std::string > deadRegionNames
 
std::vector< const G4Region * > deadRegions
 
std::vector< double > ekinMins
 
std::vector< std::string > ekinNames
 
std::vector< std::string > ekinParticles
 
std::vector< int > ekinPDG
 
std::vector< G4LogicalVolume * > ekinVolumes
 
EventActioneventAction_
 
bool initialized
 
bool killBeamPipe
 
std::vector< std::string > maxTimeNames
 
std::vector< const G4Region * > maxTimeRegions
 
double maxTrackTime
 
std::vector< double > maxTrackTimes
 
unsigned int ndeadRegions
 
unsigned int numberEkins
 
unsigned int numberPart
 
unsigned int numberTimes
 
double theCriticalDensity
 
double theCriticalEnergyForVacuum
 
G4VPhysicalVolume * tracker
 

Detailed Description

Definition at line 20 of file SteppingAction.h.

Constructor & Destructor Documentation

SteppingAction::SteppingAction ( EventAction ea,
const edm::ParameterSet ps 
)

Definition at line 16 of file SteppingAction.cc.

References deadRegionNames, ekinMins, ekinNames, ekinParticles, g, edm::ParameterSet::getParameter(), i, killBeamPipe, maxTimeNames, maxTrackTime, maxTrackTimes, MeV, ndeadRegions, numberEkins, numberPart, numberTimes, theCriticalDensity, and theCriticalEnergyForVacuum.

17  : eventAction_(e), tracker(0), calo(0), initialized(false), killBeamPipe(false)
18 {
20  (p.getParameter<double>("CriticalEnergyForVacuum")*CLHEP::MeV);
21  if(0.0 < theCriticalEnergyForVacuum) { killBeamPipe = true; }
23  (p.getParameter<double>("CriticalDensity")*CLHEP::g/CLHEP::cm3);
24  maxTrackTime = p.getParameter<double>("MaxTrackTime")*CLHEP::ns;
25  maxTrackTimes = p.getParameter<std::vector<double> >("MaxTrackTimes");
26  maxTimeNames = p.getParameter<std::vector<std::string> >("MaxTimeNames");
27  deadRegionNames = p.getParameter<std::vector<std::string> >("DeadRegions");
28  ekinMins = p.getParameter<std::vector<double> >("EkinThresholds");
29  ekinNames = p.getParameter<std::vector<std::string> >("EkinNames");
30  ekinParticles = p.getParameter<std::vector<std::string> >("EkinParticles");
31 
32  edm::LogInfo("SimG4CoreApplication") << "SteppingAction:: KillBeamPipe = "
33  << killBeamPipe << " CriticalDensity = "
34  << theCriticalDensity*CLHEP::cm3/CLHEP::g
35  << " g/cm3;"
36  << " CriticalEnergyForVacuum = "
38  << " Mev;"
39  << " MaxTrackTime = "
40  << maxTrackTime/CLHEP::ns
41  << " ns";
42 
43  numberTimes = maxTrackTimes.size();
44  if(numberTimes > 0) {
45  for (unsigned int i=0; i<numberTimes; i++) {
46  edm::LogInfo("SimG4CoreApplication") << "SteppingAction::MaxTrackTime for "
47  << maxTimeNames[i] << " is "
48  << maxTrackTimes[i] << " ns ";
49  maxTrackTimes[i] *= ns;
50  }
51  }
52 
53  ndeadRegions = deadRegionNames.size();
54  if(ndeadRegions > 0) {
55  edm::LogInfo("SimG4CoreApplication")
56  << "SteppingAction: Number of DeadRegions where all trackes are killed "
57  << ndeadRegions;
58  for(unsigned int i=0; i<ndeadRegions; ++i) {
59  edm::LogInfo("SimG4CoreApplication")
60  << "SteppingAction: DeadRegion " << i << ". " << deadRegionNames[i];
61  }
62  }
63  numberEkins = ekinNames.size();
64  numberPart = ekinParticles.size();
65  if(0 == numberPart) { numberEkins = 0; }
66 
67  if(numberEkins > 0) {
68 
69  edm::LogInfo("SimG4CoreApplication") << "SteppingAction::Kill following "
70  << numberPart
71  << " particles in " << numberEkins
72  << " volumes";
73  for (unsigned int i=0; i<numberPart; ++i) {
74  edm::LogInfo("SimG4CoreApplication") << "SteppingAction::Particle " << i
75  << " " << ekinParticles[i]
76  << " Threshold = " << ekinMins[i]
77  << " MeV";
78  ekinMins[i] *= CLHEP::MeV;
79  }
80  for (unsigned int i=0; i<numberEkins; ++i) {
81  edm::LogInfo("SimG4CoreApplication") << "SteppingAction::LogVolume[" << i
82  << "] = " << ekinNames[i];
83  }
84  }
85 }
int i
Definition: DBlmapReader.cc:9
double theCriticalDensity
G4VPhysicalVolume * calo
unsigned int numberPart
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
std::vector< std::string > deadRegionNames
unsigned int numberTimes
std::vector< double > maxTrackTimes
double theCriticalEnergyForVacuum
EventAction * eventAction_
std::vector< double > ekinMins
std::vector< std::string > maxTimeNames
std::vector< std::string > ekinNames
G4VPhysicalVolume * tracker
unsigned int ndeadRegions
std::vector< std::string > ekinParticles
SteppingAction::~SteppingAction ( )
virtual

Definition at line 87 of file SteppingAction.cc.

87 {}

Member Function Documentation

bool SteppingAction::catchLongLived ( G4Track *  theTrack,
const G4Region *  reg 
) const
private

Definition at line 166 of file SteppingAction.cc.

References archive::flag, i, maxTimeRegions, maxTrackTime, maxTrackTimes, numberTimes, and PrintKilledTrack().

Referenced by UserSteppingAction().

167 {
168  bool flag = true;
169  double tofM = maxTrackTime;
170 
171  if(numberTimes > 0) {
172  for (unsigned int i=0; i<numberTimes; ++i) {
173  if (reg == maxTimeRegions[i]) {
174  tofM = maxTrackTimes[i];
175  break;
176  }
177  }
178  }
179  if (theTrack->GetGlobalTime() > tofM) {
180  theTrack->SetTrackStatus(fStopAndKill);
181 #ifdef DebugLog
182  PrintKilledTrack(theTrack, "out of time");
183 #endif
184  flag = false;
185  }
186  return flag;
187 }
int i
Definition: DBlmapReader.cc:9
unsigned int numberTimes
std::vector< double > maxTrackTimes
void PrintKilledTrack(const G4Track *, const std::string &) const
std::vector< const G4Region * > maxTimeRegions
bool SteppingAction::initPointer ( )
private

Definition at line 232 of file SteppingAction.cc.

References calo, deadRegionNames, deadRegions, ekinMins, ekinNames, ekinParticles, ekinPDG, ekinVolumes, i, LogDebug, maxTimeNames, maxTimeRegions, MeV, ndeadRegions, numberEkins, numberPart, numberTimes, and tracker.

Referenced by UserSteppingAction().

233 {
234  const G4PhysicalVolumeStore * pvs = G4PhysicalVolumeStore::GetInstance();
235  if (pvs) {
236  std::vector<G4VPhysicalVolume*>::const_iterator pvcite;
237  for (pvcite = pvs->begin(); pvcite != pvs->end(); ++pvcite) {
238  if ((*pvcite)->GetName() == "Tracker") tracker = (*pvcite);
239  if ((*pvcite)->GetName() == "CALO") calo = (*pvcite);
240  if (tracker && calo) break;
241  }
242  if (tracker || calo) {
243  edm::LogInfo("SimG4CoreApplication") << "Pointer for Tracker " << tracker
244  << " and for Calo " << calo;
245  if (tracker) LogDebug("SimG4CoreApplication") << "Tracker vol name "
246  << tracker->GetName();
247  if (calo) LogDebug("SimG4CoreApplication") << "Calorimeter vol name "
248  << calo->GetName();
249  }
250  }
251 
252  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
253  if (numberEkins > 0) {
254  if (lvs) {
255  ekinVolumes.resize(numberEkins, 0);
256  std::vector<G4LogicalVolume*>::const_iterator lvcite;
257  for (lvcite = lvs->begin(); lvcite != lvs->end(); ++lvcite) {
258  for (unsigned int i=0; i<numberEkins; ++i) {
259  if ((*lvcite)->GetName() == (G4String)(ekinNames[i])) {
260  ekinVolumes[i] = (*lvcite);
261  break;
262  }
263  }
264  }
265  }
266  for (unsigned int i=0; i<numberEkins; ++i) {
267  edm::LogInfo("SimG4CoreApplication") << ekinVolumes[i]->GetName()
268  <<" with pointer " << ekinVolumes[i];
269  }
270  }
271 
272  if(numberPart > 0) {
273  G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
274  G4String partName;
275  ekinPDG.resize(numberPart, 0);
276  for (unsigned int i=0; i<numberPart; ++i) {
277  ekinPDG[i] =
278  theParticleTable->FindParticle(partName=ekinParticles[i])->GetPDGEncoding();
279  edm::LogInfo("SimG4CoreApplication") << "Particle " << ekinParticles[i]
280  << " with PDG code " << ekinPDG[i]
281  << " and KE cut off "
282  << ekinMins[i]/MeV << " MeV";
283  }
284  }
285 
286  const G4RegionStore * rs = G4RegionStore::GetInstance();
287  if (numberTimes > 0) {
288  maxTimeRegions.resize(numberTimes, 0);
289  std::vector<G4Region*>::const_iterator rcite;
290  for (rcite = rs->begin(); rcite != rs->end(); ++rcite) {
291  for (unsigned int i=0; i<numberTimes; ++i) {
292  if ((*rcite)->GetName() == (G4String)(maxTimeNames[i])) {
293  maxTimeRegions[i] = (*rcite);
294  break;
295  }
296  }
297  }
298  }
299  if (ndeadRegions > 0) {
300  deadRegions.resize(ndeadRegions, 0);
301  std::vector<G4Region*>::const_iterator rcite;
302  for (rcite = rs->begin(); rcite != rs->end(); ++rcite) {
303  for (unsigned int i=0; i<ndeadRegions; ++i) {
304  if ((*rcite)->GetName() == (G4String)(deadRegionNames[i])) {
305  deadRegions[i] = (*rcite);
306  break;
307  }
308  }
309  }
310  }
311  return true;
312 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
std::vector< int > ekinPDG
std::vector< const G4Region * > deadRegions
G4VPhysicalVolume * calo
unsigned int numberPart
unsigned int numberEkins
const double MeV
std::vector< std::string > deadRegionNames
unsigned int numberTimes
std::vector< G4LogicalVolume * > ekinVolumes
std::vector< double > ekinMins
std::vector< std::string > maxTimeNames
std::vector< std::string > ekinNames
G4VPhysicalVolume * tracker
unsigned int ndeadRegions
std::vector< std::string > ekinParticles
std::vector< const G4Region * > maxTimeRegions
bool SteppingAction::isThisVolume ( const G4VTouchable *  touch,
G4VPhysicalVolume *  pv 
) const
private

Definition at line 223 of file SteppingAction.cc.

References testEve_cfg::level.

Referenced by UserSteppingAction().

225 {
226  bool res = false;
227  int level = (touch->GetHistoryDepth())+1;
228  if (level >= 3) { res = (touch->GetVolume(level - 3) == pv); }
229  return res;
230 }
tuple level
Definition: testEve_cfg.py:34
bool SteppingAction::killInsideDeadRegion ( G4Track *  theTrack,
const G4Region *  reg 
) const
private

Definition at line 149 of file SteppingAction.cc.

References deadRegions, i, ndeadRegions, and PrintKilledTrack().

Referenced by UserSteppingAction().

151 {
152  bool alive = true;
153  for(unsigned int i=0; i<ndeadRegions; ++i) {
154  if(reg == deadRegions[i]) {
155  alive = false;
156  theTrack->SetTrackStatus(fStopAndKill);
157 #ifdef DebugLog
158  PrintKilledTrack(theTrack, "dead region");
159 #endif
160  break;
161  }
162  }
163  return alive;
164 }
int i
Definition: DBlmapReader.cc:9
std::vector< const G4Region * > deadRegions
void PrintKilledTrack(const G4Track *, const std::string &) const
unsigned int ndeadRegions
bool SteppingAction::killLowEnergy ( const G4Step *  aStep) const
private

Definition at line 189 of file SteppingAction.cc.

References ekinMins, ekinPDG, ekinVolumes, archive::flag, i, numberEkins, numberPart, convertSQLiteXML::ok, and PrintKilledTrack().

Referenced by UserSteppingAction().

190 {
191  bool ok = true;
192  bool flag = false;
193  G4LogicalVolume* lv =
194  aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
195  for (unsigned int i=0; i<numberEkins; ++i) {
196  if (lv == ekinVolumes[i]) {
197  flag = true;
198  break;
199  }
200  }
201  if (flag) {
202  G4Track * track = aStep->GetTrack();
203  double ekin = track->GetKineticEnergy();
204  double ekinM = 0;
205  int pCode = track->GetDefinition()->GetPDGEncoding();
206  for (unsigned int i=0; i<numberPart; ++i) {
207  if (pCode == ekinPDG[i]) {
208  ekinM = ekinMins[i];
209  break;
210  }
211  }
212  if (ekin <= ekinM) {
213  track->SetTrackStatus(fStopAndKill);
214 #ifdef DebugLog
215  PrintKilledTrack(track, "low-energy");
216 #endif
217  ok = false;
218  }
219  }
220  return ok;
221 }
int i
Definition: DBlmapReader.cc:9
std::vector< int > ekinPDG
unsigned int numberPart
unsigned int numberEkins
std::vector< G4LogicalVolume * > ekinVolumes
std::vector< double > ekinMins
void PrintKilledTrack(const G4Track *, const std::string &) const
void SteppingAction::PrintKilledTrack ( const G4Track *  aTrack,
const std::string &  typ 
) const
private

Definition at line 314 of file SteppingAction.cc.

References MeV, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by catchLongLived(), killInsideDeadRegion(), killLowEnergy(), and UserSteppingAction().

316 {
317  std::string vname = "";
318  std::string rname = "";
319  G4VPhysicalVolume* pv = aTrack->GetNextVolume();
320  if(pv) {
321  vname = pv->GetLogicalVolume()->GetName();
322  rname = pv->GetLogicalVolume()->GetRegion()->GetName();
323  }
324 
325  edm::LogInfo("SimG4CoreApplication")
326  << "Track #" << aTrack->GetTrackID()
327  << " " << aTrack->GetDefinition()->GetParticleName()
328  << " E(MeV)= " << aTrack->GetKineticEnergy()/MeV
329  << " is killed due to " << typ
330  << " inside LV: " << vname << " (" << rname
331  << ") at " << aTrack->GetPosition();
332 }
const double MeV
void SteppingAction::UserSteppingAction ( const G4Step *  aStep)
virtual

Definition at line 89 of file SteppingAction.cc.

References EventAction::addTkCaloStateInfo(), calo, catchLongLived(), eventAction_, initialized, initPointer(), isThisVolume(), killBeamPipe, killInsideDeadRegion(), killLowEnergy(), m_g4StepSignal, ndeadRegions, numberEkins, convertSQLiteXML::ok, AlCaHLTBitMon_ParallelJobs::p, PrintKilledTrack(), theCriticalDensity, theCriticalEnergyForVacuum, tracker, x, detailsBasic3DVector::y, and detailsBasic3DVector::z.

90 {
91  if (!initialized) { initialized = initPointer(); }
92  m_g4StepSignal(aStep);
93 
94  G4Track * theTrack = aStep->GetTrack();
95  bool ok = (theTrack->GetTrackStatus() == fAlive);
96  G4StepPoint* postStep = aStep->GetPostStepPoint();
97  if(ok && postStep->GetPhysicalVolume() != 0) {
98 
99  G4StepPoint* preStep = aStep->GetPreStepPoint();
100  const G4Region* theRegion =
101  preStep->GetPhysicalVolume()->GetLogicalVolume()->GetRegion();
102 
103  // kill in dead regions
104  if(ok && 0 < ndeadRegions) { ok = killInsideDeadRegion(theTrack, theRegion); }
105 
106  // kill out of time
107  if(ok) { ok = catchLongLived(theTrack, theRegion); }
108 
109  // kill low-energy in volumes on demand
110  if(ok && numberEkins > 0) { ok = killLowEnergy(aStep); }
111 
112  // kill low-energy in vacuum
113  G4double kinEnergy = theTrack->GetKineticEnergy();
114  if(ok && killBeamPipe && kinEnergy < theCriticalEnergyForVacuum
115  && theTrack->GetDefinition()->GetPDGCharge() != 0.0 && kinEnergy > 0.0
116  && theTrack->GetNextVolume()->GetLogicalVolume()->GetMaterial()->GetDensity()
117  <= theCriticalDensity) {
118  theTrack->SetTrackStatus(fStopAndKill);
119 #ifdef DebugLog
120  PrintKilledTrack(theTrack, "LE in vacuum");
121 #endif
122  ok = false;
123  }
124 
125  // check transition tracker/calo
126  if(ok) {
127 
128  if(isThisVolume(preStep->GetTouchable(),tracker) &&
129  isThisVolume(postStep->GetTouchable(),calo)) {
130 
131  math::XYZVectorD pos((preStep->GetPosition()).x(),
132  (preStep->GetPosition()).y(),
133  (preStep->GetPosition()).z());
134 
135  math::XYZTLorentzVectorD mom((preStep->GetMomentum()).x(),
136  (preStep->GetMomentum()).y(),
137  (preStep->GetMomentum()).z(),
138  preStep->GetTotalEnergy());
139 
140  uint32_t id = theTrack->GetTrackID();
141 
142  std::pair<math::XYZVectorD,math::XYZTLorentzVectorD> p(pos,mom);
144  }
145  }
146  }
147 }
bool killInsideDeadRegion(G4Track *theTrack, const G4Region *reg) const
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
bool killLowEnergy(const G4Step *aStep) const
double theCriticalDensity
SimActivityRegistry::G4StepSignal m_g4StepSignal
G4VPhysicalVolume * calo
unsigned int numberEkins
float float float z
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
double theCriticalEnergyForVacuum
EventAction * eventAction_
bool isThisVolume(const G4VTouchable *touch, G4VPhysicalVolume *pv) const
void PrintKilledTrack(const G4Track *, const std::string &) const
void addTkCaloStateInfo(uint32_t t, const std::pair< math::XYZVectorD, math::XYZTLorentzVectorD > &p)
Definition: EventAction.cc:90
G4VPhysicalVolume * tracker
Definition: DDAxes.h:10
unsigned int ndeadRegions
bool catchLongLived(G4Track *theTrack, const G4Region *reg) const

Member Data Documentation

G4VPhysicalVolume * SteppingAction::calo
private

Definition at line 44 of file SteppingAction.h.

Referenced by initPointer(), and UserSteppingAction().

std::vector<std::string> SteppingAction::deadRegionNames
private

Definition at line 50 of file SteppingAction.h.

Referenced by initPointer(), and SteppingAction().

std::vector<const G4Region*> SteppingAction::deadRegions
private

Definition at line 52 of file SteppingAction.h.

Referenced by initPointer(), and killInsideDeadRegion().

std::vector<double> SteppingAction::ekinMins
private

Definition at line 48 of file SteppingAction.h.

Referenced by initPointer(), killLowEnergy(), and SteppingAction().

std::vector<std::string> SteppingAction::ekinNames
private

Definition at line 49 of file SteppingAction.h.

Referenced by initPointer(), and SteppingAction().

std::vector<std::string> SteppingAction::ekinParticles
private

Definition at line 49 of file SteppingAction.h.

Referenced by initPointer(), and SteppingAction().

std::vector<int> SteppingAction::ekinPDG
private

Definition at line 54 of file SteppingAction.h.

Referenced by initPointer(), and killLowEnergy().

std::vector<G4LogicalVolume*> SteppingAction::ekinVolumes
private

Definition at line 53 of file SteppingAction.h.

Referenced by initPointer(), and killLowEnergy().

EventAction* SteppingAction::eventAction_
private

Definition at line 43 of file SteppingAction.h.

Referenced by UserSteppingAction().

bool SteppingAction::initialized
private

Definition at line 60 of file SteppingAction.h.

Referenced by UserSteppingAction().

bool SteppingAction::killBeamPipe
private

Definition at line 61 of file SteppingAction.h.

Referenced by SteppingAction(), and UserSteppingAction().

SimActivityRegistry::G4StepSignal SteppingAction::m_g4StepSignal
std::vector<std::string> SteppingAction::maxTimeNames
private

Definition at line 49 of file SteppingAction.h.

Referenced by initPointer(), and SteppingAction().

std::vector<const G4Region*> SteppingAction::maxTimeRegions
private

Definition at line 51 of file SteppingAction.h.

Referenced by catchLongLived(), and initPointer().

double SteppingAction::maxTrackTime
private

Definition at line 47 of file SteppingAction.h.

Referenced by catchLongLived(), and SteppingAction().

std::vector<double> SteppingAction::maxTrackTimes
private

Definition at line 48 of file SteppingAction.h.

Referenced by catchLongLived(), and SteppingAction().

unsigned int SteppingAction::ndeadRegions
private
unsigned int SteppingAction::numberEkins
private

Definition at line 56 of file SteppingAction.h.

Referenced by initPointer(), killLowEnergy(), SteppingAction(), and UserSteppingAction().

unsigned int SteppingAction::numberPart
private

Definition at line 57 of file SteppingAction.h.

Referenced by initPointer(), killLowEnergy(), and SteppingAction().

unsigned int SteppingAction::numberTimes
private

Definition at line 55 of file SteppingAction.h.

Referenced by catchLongLived(), initPointer(), and SteppingAction().

double SteppingAction::theCriticalDensity
private

Definition at line 46 of file SteppingAction.h.

Referenced by SteppingAction(), and UserSteppingAction().

double SteppingAction::theCriticalEnergyForVacuum
private

Definition at line 45 of file SteppingAction.h.

Referenced by SteppingAction(), and UserSteppingAction().

G4VPhysicalVolume* SteppingAction::tracker
private

Definition at line 44 of file SteppingAction.h.

Referenced by initPointer(), and UserSteppingAction().