CMS 3D CMS Logo

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

#include <SimG4CMS/Forward/interface/TotemSD.h>

Inheritance diagram for TotemSD:
SensitiveTkDetector Observer< const BeginOfEvent * > Observer< const EndOfEvent * > SensitiveDetector

Public Member Functions

virtual void clear ()
 
virtual void DrawAll ()
 
virtual void EndOfEvent (G4HCofThisEvent *eventHC)
 
void fillHits (edm::PSimHitContainer &, std::string use)
 
virtual void Initialize (G4HCofThisEvent *HCE)
 
virtual void PrintAll ()
 
virtual bool ProcessHits (G4Step *, G4TouchableHistory *)
 
virtual uint32_t setDetUnitId (G4Step *)
 
 TotemSD (std::string, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
virtual ~TotemSD ()
 
- Public Member Functions inherited from SensitiveTkDetector
 SensitiveTkDetector (std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
 
- Public Member Functions inherited from SensitiveDetector
virtual void AssignSD (std::string &vname)
 
Local3DPoint ConvertToLocal3DPoint (G4ThreeVector point)
 
Local3DPoint FinalStepPosition (G4Step *s, coordinates)
 
virtual std::vector< std::string > getNames ()
 
Local3DPoint InitialStepPosition (G4Step *s, coordinates)
 
std::string nameOfSD ()
 
void NaNTrap (G4Step *step)
 
void Register ()
 
 SensitiveDetector (std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &, edm::ParameterSet const &p)
 
virtual ~SensitiveDetector ()
 
- Public Member Functions inherited from Observer< const BeginOfEvent * >
 Observer ()
 
void slotForUpdate (const BeginOfEvent *iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfEvent * >
 Observer ()
 
void slotForUpdate (const EndOfEvent *iT)
 
virtual ~Observer ()
 

Private Member Functions

virtual void clearHits ()
 
void CreateNewHit ()
 
void CreateNewHitEvo ()
 
void GetStepInfo (G4Step *aStep)
 
bool HitExists ()
 
G4ThreeVector PosizioEvo (G4ThreeVector, double, double, double, double, int &)
 
void ResetForNewPrimary ()
 
G4ThreeVector SetToLocal (G4ThreeVector globalPoint)
 
void StoreHit (TotemG4Hit *)
 
void Summarize ()
 
void update (const BeginOfEvent *)
 This routine will be called when the appropriate signal arrives. More...
 
void update (const ::EndOfEvent *)
 
void UpdateHit ()
 

Private Attributes

TotemG4HitcurrentHit
 
G4VPhysicalVolume * currentPV
 
float edeposit
 
float Eloss
 
G4ThreeVector entrancePoint
 
int eventno
 
G4int hcID
 
G4ThreeVector hitPoint
 
float incidentEnergy
 
std::string name
 
TotemVDetectorOrganizationnumberingScheme
 
float Pabs
 
int ParentId
 
short ParticleType
 
float PhiAtEntry
 
G4ThreeVector Posizio
 
G4StepPoint * postStepPoint
 
G4StepPoint * preStepPoint
 
uint32_t previousUnitID
 
int primaryID
 
G4int primID
 
TrackingSlaveSDslave
 
TotemG4HitCollectiontheHC
 
const SimTrackManagertheManager
 
float ThetaAtEntry
 
G4Track * theTrack
 
float Tof
 
int tsID
 
double tSlice
 
int tSliceID
 
uint32_t unitID
 
float Vx
 
float Vy
 
float Vz
 

Additional Inherited Members

- Public Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 
- Protected Member Functions inherited from Observer< const EndOfEvent * >
virtual void update (const EndOfEvent *)=0
 This routine will be called when the appropriate signal arrives. More...
 

Detailed Description

Description: Stores hits of Totem in appropriate container

Usage: Used in sensitive detector builder

Definition at line 43 of file TotemSD.h.

Constructor & Destructor Documentation

TotemSD::TotemSD ( std::string  name,
const DDCompactView cpv,
SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)

Definition at line 45 of file TotemSD.cc.

References SensitiveDetector::AssignSD(), ecaldqm::collectionName, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), LogDebug, SensitiveDetectorCatalog::logicalNames(), numberingScheme, SensitiveDetector::Register(), and slave.

47  :
49  hcID(-1), theHC(0), theManager(manager), currentHit(0), theTrack(0),
51  postStepPoint(0), eventno(0){
52 
53  //Add Totem Sentitive Detector Names
54  collectionName.insert(name);
55 
56  //Parameters
57  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("TotemSD");
58  int verbn = m_p.getUntrackedParameter<int>("Verbosity");
59 
60  SetVerboseLevel(verbn);
61  LogDebug("ForwardSim")
62  << "*******************************************************\n"
63  << "* *\n"
64  << "* Constructing a TotemSD with name " << name << "\n"
65  << "* *\n"
66  << "*******************************************************";
67 
68  slave = new TrackingSlaveSD(name);
69 
70  //
71  // Now attach the right detectors (LogicalVolumes) to me
72  //
73  std::vector<std::string> lvNames = clg.logicalNames(name);
74  this->Register();
75  for (std::vector<std::string>::iterator it=lvNames.begin();
76  it !=lvNames.end(); it++) {
77  this->AssignSD(*it);
78  edm::LogInfo("ForwardSim") << "TotemSD : Assigns SD to LV " << (*it);
79  }
80 
81  if (name == "TotemHitsT1") {
83  } else if (name == "TotemHitsT2Si") {
85  } else if (name == "TotemHitsT2Gem") {
87  } else if (name == "TotemHitsRP") {
89  } else {
90  edm::LogWarning("ForwardSim") << "TotemSD: ReadoutName not supported\n";
91  }
92 
93  edm::LogInfo("ForwardSim") << "TotemSD: Instantiation completed";
94 }
#define LogDebug(id)
T getUntrackedParameter(std::string const &, T const &) const
uint32_t unitID
Definition: TotemSD.h:106
std::vector< std::string > logicalNames(std::string &readoutName)
TotemG4HitCollection * theHC
Definition: TotemSD.h:99
uint32_t previousUnitID
Definition: TotemSD.h:106
TotemVDetectorOrganization * numberingScheme
Definition: TotemSD.h:86
G4int hcID
Definition: TotemSD.h:98
G4StepPoint * preStepPoint
Definition: TotemSD.h:110
SensitiveTkDetector(std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
std::string const collectionName[nCollections]
Definition: Collections.h:39
G4StepPoint * postStepPoint
Definition: TotemSD.h:111
TotemG4Hit * currentHit
Definition: TotemSD.h:103
G4VPhysicalVolume * currentPV
Definition: TotemSD.h:105
virtual void AssignSD(std::string &vname)
TrackingSlaveSD * slave
Definition: TotemSD.h:85
G4Track * theTrack
Definition: TotemSD.h:104
int eventno
Definition: TotemSD.h:127
const SimTrackManager * theManager
Definition: TotemSD.h:100
std::string name
Definition: TotemSD.h:97
TotemSD::~TotemSD ( )
virtual

Definition at line 96 of file TotemSD.cc.

References numberingScheme, and slave.

96  {
97  if (slave) delete slave;
98  if (numberingScheme) delete numberingScheme;
99 }
TotemVDetectorOrganization * numberingScheme
Definition: TotemSD.h:86
TrackingSlaveSD * slave
Definition: TotemSD.h:85

Member Function Documentation

void TotemSD::clear ( void  )
virtual

Definition at line 165 of file TotemSD.cc.

165  {
166 }
void TotemSD::clearHits ( )
privatevirtual

Implements SensitiveDetector.

Definition at line 190 of file TotemSD.cc.

References TrackingSlaveSD::Initialize(), and slave.

Referenced by update().

190  {
191  slave->Initialize();
192 }
virtual void Initialize()
TrackingSlaveSD * slave
Definition: TotemSD.h:85
void TotemSD::CreateNewHit ( )
private

Definition at line 285 of file TotemSD.cc.

References gather_cfg::cout, currentHit, currentPV, Eloss, incidentEnergy, LogDebug, NULL, Pabs, ParentId, ParticleType, PhiAtEntry, Posizio, primaryID, TotemG4Hit::setEnergyLoss(), TotemG4Hit::setEntry(), TotemG4Hit::setIncidentEnergy(), TotemG4Hit::setPabs(), TotemG4Hit::setParentId(), TotemG4Hit::setParticleType(), TotemG4Hit::setPhiAtEntry(), TotemG4Hit::setThetaAtEntry(), TotemG4Hit::setTimeSlice(), TotemG4Hit::setTof(), TotemG4Hit::setTrackID(), TotemG4Hit::setUnitID(), TotemG4Hit::setVx(), TotemG4Hit::setVy(), TotemG4Hit::setVz(), StoreHit(), ThetaAtEntry, theTrack, Tof, tSlice, tSliceID, unitID, UpdateHit(), Vx, Vy, and Vz.

Referenced by ProcessHits().

285  {
286 
287 #ifdef debug
288  LogDebug("ForwardSim") << "TotemSD CreateNewHit for"
289  << " PV " << currentPV->GetName()
290  << " PVid = " << currentPV->GetCopyNo()
291  << " MVid = " << currentPV->GetMother()->GetCopyNo()
292  << " Unit " << unitID << "\n"
293  << " primary " << primaryID
294  << " time slice " << tSliceID
295  << " For Track " << theTrack->GetTrackID()
296  << " which is a "
297  << theTrack->GetDefinition()->GetParticleName();
298 
299  if (theTrack->GetTrackID()==1) {
300  LogDebug("ForwardSim") << " of energy " << theTrack->GetTotalEnergy();
301  } else {
302  LogDebug("ForwardSim") << " daughter of part. " << theTrack->GetParentID();
303  }
304 
305  cout << " and created by " ;
306  if (theTrack->GetCreatorProcess()!=NULL)
307  LogDebug("ForwardSim") << theTrack->GetCreatorProcess()->GetProcessName() ;
308  else
309  LogDebug("ForwardSim") << "NO process";
310 #endif
311 
312 
313  currentHit = new TotemG4Hit;
318 
325 
326  currentHit->setEntry(Posizio.x(),Posizio.y(),Posizio.z());
327 
329  currentHit->setVx(Vx);
330  currentHit->setVy(Vy);
331  currentHit->setVz(Vz);
332 
333  UpdateHit();
334 
336 }
#define LogDebug(id)
void setVz(float p)
Definition: TotemG4Hit.cc:184
int ParentId
Definition: TotemSD.h:124
void setPabs(float e)
Definition: TotemG4Hit.cc:154
uint32_t unitID
Definition: TotemSD.h:106
float Pabs
Definition: TotemSD.h:116
G4ThreeVector Posizio
Definition: TotemSD.h:115
float Tof
Definition: TotemSD.h:117
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:156
void setEntry(double x, double y, double z)
Definition: TotemG4Hit.h:58
double tSlice
Definition: TotemSD.h:108
int primaryID
Definition: TotemSD.h:107
float incidentEnergy
Definition: TotemSD.h:94
#define NULL
Definition: scimark2.h:8
int tSliceID
Definition: TotemSD.h:107
void setVy(float p)
Definition: TotemG4Hit.cc:181
void setTimeSlice(double d)
Definition: TotemG4Hit.cc:142
float Eloss
Definition: TotemSD.h:118
float Vx
Definition: TotemSD.h:125
void setVx(float p)
Definition: TotemG4Hit.cc:178
void setParentId(int p)
Definition: TotemG4Hit.cc:175
TotemG4Hit * currentHit
Definition: TotemSD.h:103
G4VPhysicalVolume * currentPV
Definition: TotemSD.h:105
void setIncidentEnergy(double e)
Definition: TotemG4Hit.cc:133
short ParticleType
Definition: TotemSD.h:119
void setThetaAtEntry(float t)
Definition: TotemG4Hit.cc:162
float PhiAtEntry
Definition: TotemSD.h:122
void setUnitID(uint32_t i)
Definition: TotemG4Hit.cc:139
void setTrackID(int i)
Definition: TotemG4Hit.cc:136
void setPhiAtEntry(float f)
Definition: TotemG4Hit.cc:163
void setParticleType(short i)
Definition: TotemG4Hit.cc:157
void setTof(float e)
Definition: TotemG4Hit.cc:155
float Vz
Definition: TotemSD.h:125
void StoreHit(TotemG4Hit *)
Definition: TotemSD.cc:505
G4Track * theTrack
Definition: TotemSD.h:104
void UpdateHit()
Definition: TotemSD.cc:473
tuple cout
Definition: gather_cfg.py:121
float Vy
Definition: TotemSD.h:125
float ThetaAtEntry
Definition: TotemSD.h:121
void TotemSD::CreateNewHitEvo ( )
private

Definition at line 338 of file TotemSD.cc.

References currentHit, Eloss, incidentEnergy, Pabs, ParentId, ParticleType, PhiAtEntry, Posizio, PosizioEvo(), primaryID, TotemG4Hit::setEnergyLoss(), TotemG4Hit::setEntry(), TotemG4Hit::setIncidentEnergy(), TotemG4Hit::setPabs(), TotemG4Hit::setParentId(), TotemG4Hit::setParticleType(), TotemG4Hit::setPhiAtEntry(), TotemG4Hit::setThetaAtEntry(), TotemG4Hit::setTimeSlice(), TotemG4Hit::setTof(), TotemG4Hit::setTrackID(), TotemG4Hit::setUnitID(), TotemG4Hit::setVx(), TotemG4Hit::setVy(), TotemG4Hit::setVz(), StoreHit(), ThetaAtEntry, Tof, tSlice, unitID, UpdateHit(), Vx, Vy, and Vz.

Referenced by ProcessHits().

338  {
339 
340 // LogDebug("ForwardSim") << "INSIDE CREATE NEW HIT EVO ";
341 
342  currentHit = new TotemG4Hit;
347 
354 
355  // LogDebug("ForwardSim") << Posizio.x() << " " << Posizio.y() << " " << Posizio.z();
356 
358  currentHit->setVx(Vx);
359  currentHit->setVy(Vy);
360  currentHit->setVz(Vz);
361 
362  G4ThreeVector _PosizioEvo;
363  int flagAcc=0;
364  _PosizioEvo=PosizioEvo(Posizio,Vx,Vy,Vz,Pabs,flagAcc);
365 
366  if(flagAcc==1){
367  currentHit->setEntry(_PosizioEvo.x(),_PosizioEvo.y(),_PosizioEvo.z());
368 
369  // if(flagAcc==1)
370  UpdateHit();
371 
373  }
374  // LogDebug("ForwardSim") << "STORED HIT IN: " << unitID;
375 }
void setVz(float p)
Definition: TotemG4Hit.cc:184
int ParentId
Definition: TotemSD.h:124
void setPabs(float e)
Definition: TotemG4Hit.cc:154
uint32_t unitID
Definition: TotemSD.h:106
float Pabs
Definition: TotemSD.h:116
G4ThreeVector Posizio
Definition: TotemSD.h:115
float Tof
Definition: TotemSD.h:117
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:156
void setEntry(double x, double y, double z)
Definition: TotemG4Hit.h:58
double tSlice
Definition: TotemSD.h:108
int primaryID
Definition: TotemSD.h:107
float incidentEnergy
Definition: TotemSD.h:94
void setVy(float p)
Definition: TotemG4Hit.cc:181
void setTimeSlice(double d)
Definition: TotemG4Hit.cc:142
float Eloss
Definition: TotemSD.h:118
float Vx
Definition: TotemSD.h:125
void setVx(float p)
Definition: TotemG4Hit.cc:178
void setParentId(int p)
Definition: TotemG4Hit.cc:175
TotemG4Hit * currentHit
Definition: TotemSD.h:103
void setIncidentEnergy(double e)
Definition: TotemG4Hit.cc:133
short ParticleType
Definition: TotemSD.h:119
void setThetaAtEntry(float t)
Definition: TotemG4Hit.cc:162
float PhiAtEntry
Definition: TotemSD.h:122
void setUnitID(uint32_t i)
Definition: TotemG4Hit.cc:139
void setTrackID(int i)
Definition: TotemG4Hit.cc:136
void setPhiAtEntry(float f)
Definition: TotemG4Hit.cc:163
void setParticleType(short i)
Definition: TotemG4Hit.cc:157
void setTof(float e)
Definition: TotemG4Hit.cc:155
float Vz
Definition: TotemSD.h:125
void StoreHit(TotemG4Hit *)
Definition: TotemSD.cc:505
G4ThreeVector PosizioEvo(G4ThreeVector, double, double, double, double, int &)
Definition: TotemSD.cc:377
void UpdateHit()
Definition: TotemSD.cc:473
float Vy
Definition: TotemSD.h:125
float ThetaAtEntry
Definition: TotemSD.h:121
void TotemSD::DrawAll ( )
virtual

Definition at line 168 of file TotemSD.cc.

168  {
169 }
void TotemSD::EndOfEvent ( G4HCofThisEvent *  eventHC)
virtual

Reimplemented from SensitiveDetector.

Definition at line 138 of file TotemSD.cc.

References TotemG4Hit::getEnergyLoss(), TotemG4Hit::getEntry(), TotemG4Hit::getPabs(), TotemG4Hit::getParticleType(), TotemG4Hit::getPhiAtEntry(), TotemG4Hit::getThetaAtEntry(), TotemG4Hit::getTof(), TotemG4Hit::getTrackID(), TotemG4Hit::getUnitID(), j, LogDebug, TrackingSlaveSD::processHits(), slave, Summarize(), and theHC.

138  {
139 
140  // here we loop over transient hits and make them persistent
141  for (int j=0; j<theHC->entries() && j<15000; j++) {
142  TotemG4Hit* aHit = (*theHC)[j];
143 #ifdef ddebug
144  LogDebug("ForwardSim") << "HIT NUMERO " << j << "unit ID = "
145  << aHit->getUnitID() << "\n"
146  << " " << "enrty z "
147  << aHit->getEntry().z() << "\n"
148  << " " << "theta "
149  << aHit->getThetaAtEntry() << "\n";
150 #endif
151  Local3DPoint theExitPoint(0,0,0);
152  Local3DPoint Entrata(aHit->getEntry().x(),
153  aHit->getEntry().y(),
154  aHit->getEntry().z());
155  slave->processHits(PSimHit(Entrata,theExitPoint,
156  aHit->getPabs(), aHit->getTof(),
157  aHit->getEnergyLoss(), aHit->getParticleType(),
158  aHit->getUnitID(), aHit->getTrackID(),
159  aHit->getThetaAtEntry(),aHit->getPhiAtEntry()));
160 
161  }
162  Summarize();
163 }
#define LogDebug(id)
TotemG4HitCollection * theHC
Definition: TotemSD.h:99
float getTof() const
Definition: TotemG4Hit.cc:150
float getEnergyLoss() const
Definition: TotemG4Hit.cc:151
math::XYZPoint getEntry() const
Definition: TotemG4Hit.cc:124
float getPabs() const
Definition: TotemG4Hit.cc:149
float getPhiAtEntry() const
Definition: TotemG4Hit.cc:160
int j
Definition: DBlmapReader.cc:9
int getTrackID() const
Definition: TotemG4Hit.cc:135
void Summarize()
Definition: TotemSD.cc:522
uint32_t getUnitID() const
Definition: TotemG4Hit.cc:138
TrackingSlaveSD * slave
Definition: TotemSD.h:85
int getParticleType() const
Definition: TotemG4Hit.cc:152
virtual bool processHits(const PSimHit &)
float getThetaAtEntry() const
Definition: TotemG4Hit.cc:159
void TotemSD::fillHits ( edm::PSimHitContainer c,
std::string  use 
)
virtual

Implements SensitiveTkDetector.

Definition at line 176 of file TotemSD.cc.

References TrackingSlaveSD::hits(), n, TrackingSlaveSD::name(), and slave.

176  {
177  if (slave->name() == n) c=slave->hits();
178 }
std::string name() const
std::vector< PSimHit > & hits()
TrackingSlaveSD * slave
Definition: TotemSD.h:85
void TotemSD::GetStepInfo ( G4Step *  aStep)
private

Definition at line 202 of file TotemSD.cc.

References currentPV, edeposit, Eloss, hitPoint, LogDebug, name, Pabs, ParentId, ParticleType, PhiAtEntry, Posizio, postStepPoint, preStepPoint, primaryID, setDetUnitId(), ThetaAtEntry, theTrack, Tof, tSlice, tSliceID, unitID, Vx, Vy, and Vz.

Referenced by ProcessHits().

202  {
203 
204  preStepPoint = aStep->GetPreStepPoint();
205  postStepPoint= aStep->GetPostStepPoint();
206  theTrack = aStep->GetTrack();
207  //Local3DPoint theEntryPoint = SensitiveDetector::InitialStepPosition(aStep,LocalCoordinates);
208  //Local3DPoint theExitPoint = SensitiveDetector::FinalStepPosition(aStep,LocalCoordinates);
209  hitPoint = preStepPoint->GetPosition();
210  currentPV = preStepPoint->GetPhysicalVolume();
211 
212  // double weight = 1;
213  G4String name = currentPV->GetName();
214  name.assign(name,0,4);
215  G4String particleType = theTrack->GetDefinition()->GetParticleName();
216  edeposit = aStep->GetTotalEnergyDeposit();
217 
218  tSlice = (postStepPoint->GetGlobalTime() )/nanosecond;
219  tSliceID = (int) tSlice;
220  unitID = setDetUnitId(aStep);
221 #ifdef debug
222  LogDebug("ForwardSim") << "UNITa " << unitID;
223 #endif
224  primaryID = theTrack->GetTrackID();
225 
226 
227  Posizio = hitPoint;
228  Pabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
229  Tof = aStep->GetPostStepPoint()->GetGlobalTime()/nanosecond;
230 
231  Eloss = aStep->GetTotalEnergyDeposit()/GeV;
232  ParticleType = theTrack->GetDefinition()->GetPDGEncoding();
233 
234  ThetaAtEntry = aStep->GetPreStepPoint()->GetPosition().theta()/deg;
235  PhiAtEntry = aStep->GetPreStepPoint()->GetPosition().phi()/deg;
236 
237  ParentId = theTrack->GetParentID();
238  Vx = theTrack->GetVertexPosition().x();
239  Vy = theTrack->GetVertexPosition().y();
240  Vz = theTrack->GetVertexPosition().z();
241 }
#define LogDebug(id)
int ParentId
Definition: TotemSD.h:124
uint32_t unitID
Definition: TotemSD.h:106
float Pabs
Definition: TotemSD.h:116
G4ThreeVector Posizio
Definition: TotemSD.h:115
float Tof
Definition: TotemSD.h:117
double tSlice
Definition: TotemSD.h:108
int primaryID
Definition: TotemSD.h:107
virtual uint32_t setDetUnitId(G4Step *)
Definition: TotemSD.cc:120
int tSliceID
Definition: TotemSD.h:107
G4StepPoint * preStepPoint
Definition: TotemSD.h:110
float Eloss
Definition: TotemSD.h:118
float edeposit
Definition: TotemSD.h:112
float Vx
Definition: TotemSD.h:125
G4ThreeVector hitPoint
Definition: TotemSD.h:113
G4StepPoint * postStepPoint
Definition: TotemSD.h:111
G4VPhysicalVolume * currentPV
Definition: TotemSD.h:105
short ParticleType
Definition: TotemSD.h:119
float PhiAtEntry
Definition: TotemSD.h:122
float Vz
Definition: TotemSD.h:125
G4Track * theTrack
Definition: TotemSD.h:104
float Vy
Definition: TotemSD.h:125
float ThetaAtEntry
Definition: TotemSD.h:121
std::string name
Definition: TotemSD.h:97
bool TotemSD::HitExists ( )
private

Definition at line 243 of file TotemSD.cc.

References currentHit, newFWLiteAna::found, TotemG4Hit::getTimeSliceID(), TotemG4Hit::getTrackID(), TotemG4Hit::getUnitID(), j, previousUnitID, primaryID, primID, ResetForNewPrimary(), theHC, tsID, tSliceID, unitID, and UpdateHit().

Referenced by ProcessHits().

243  {
244 
245  if (primaryID<1) {
246  edm::LogWarning("ForwardSim") << "***** TotemSD error: primaryID = "
247  << primaryID
248  << " maybe detector name changed";
249  }
250 
251  // Update if in the same detector, time-slice and for same track
252  // if (primaryID == primID && tSliceID == tsID && unitID==previousUnitID) {
253  if (tSliceID == tsID && unitID==previousUnitID) {
254  UpdateHit();
255  return true;
256  }
257 
258  // Reset entry point for new primary
259  if (primaryID != primID)
261 
262  //look in the HitContainer whether a hit with the same primID, unitID,
263  //tSliceID already exists:
264 
265  bool found = false;
266 
267  for (int j=0; j<theHC->entries()&&!found; j++) {
268  TotemG4Hit* aPreviousHit = (*theHC)[j];
269  if (aPreviousHit->getTrackID() == primaryID &&
270  aPreviousHit->getTimeSliceID() == tSliceID &&
271  aPreviousHit->getUnitID() == unitID ) {
272  currentHit = aPreviousHit;
273  found = true;
274  }
275  }
276 
277  if (found) {
278  UpdateHit();
279  return true;
280  } else {
281  return false;
282  }
283 }
uint32_t unitID
Definition: TotemSD.h:106
TotemG4HitCollection * theHC
Definition: TotemSD.h:99
int primaryID
Definition: TotemSD.h:107
uint32_t previousUnitID
Definition: TotemSD.h:106
int tSliceID
Definition: TotemSD.h:107
void ResetForNewPrimary()
Definition: TotemSD.cc:516
int j
Definition: DBlmapReader.cc:9
int getTrackID() const
Definition: TotemG4Hit.cc:135
int getTimeSliceID() const
Definition: TotemG4Hit.cc:143
TotemG4Hit * currentHit
Definition: TotemSD.h:103
uint32_t getUnitID() const
Definition: TotemG4Hit.cc:138
G4int primID
Definition: TotemSD.h:95
void UpdateHit()
Definition: TotemSD.cc:473
int tsID
Definition: TotemSD.h:102
void TotemSD::Initialize ( G4HCofThisEvent *  HCE)
virtual

Reimplemented from SensitiveDetector.

Definition at line 125 of file TotemSD.cc.

References ecaldqm::collectionName, hcID, LogDebug, name, primID, theHC, and tsID.

125  {
126 
127  LogDebug("ForwardSim") << "TotemSD : Initialize called for " << name;
128 
130  if (hcID<0)
131  hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
132  HCE->AddHitsCollection(hcID, theHC);
133 
134  tsID = -2;
135  primID = -2;
136 }
#define LogDebug(id)
TotemG4HitCollection * theHC
Definition: TotemSD.h:99
G4int hcID
Definition: TotemSD.h:98
std::string const collectionName[nCollections]
Definition: Collections.h:39
G4int primID
Definition: TotemSD.h:95
int tsID
Definition: TotemSD.h:102
G4THitsCollection< TotemG4Hit > TotemG4HitCollection
std::string name
Definition: TotemSD.h:97
G4ThreeVector TotemSD::PosizioEvo ( G4ThreeVector  Pos,
double  vx,
double  vy,
double  vz,
double  pabs,
int &  accettanza 
)
private

Definition at line 377 of file TotemSD.cc.

References j, and mathSSE::sqrt().

Referenced by CreateNewHitEvo().

378  {
379  accettanza=0;
380  //Pos.xyz() in mm
381  G4ThreeVector PosEvo;
382  double ThetaX=atan((Pos.x()-vx)/(Pos.z()-vz));
383  double ThetaY=atan((Pos.y()-vy)/(Pos.z()-vz));
384  double X_at_0 =(vx-((Pos.x()-vx)/(Pos.z()-vz))*vz)/1000.;
385  double Y_at_0 =(vy-((Pos.y()-vy)/(Pos.z()-vz))*vz)/1000.;
386 
387  // double temp_evo_X;
388  // double temp_evo_Y;
389  // double temp_evo_Z;
390  // temp_evo_Z = fabs(Pos.z())/Pos.z()*220000.;
391 
392  //csi=-dp/d
393  double csi = fabs((7000.-pabs)/7000.);
394 
395  // all in m
396  const int no_rp=4;
397  double x_par[no_rp+1];
398  double y_par[no_rp+1];
399  //rp z position
400  double rp[no_rp]={141.,149.,198.,220.};
401  //{lx0,mlx} for each rp; Lx=lx0+mlx*csi
402  double leffx[][2]={{122.5429,-46.9312},{125.4194,-49.1849},{152.6,-81.157},{98.8914,-131.8390}};
403  //{ly0,mly} for each rp; Ly=ly0+mly*csi
404  double leffy[][2]={{124.2314,-55.4852},{127.7825,-57.4503},{179.455,-76.274},{273.0931,-40.4626}};
405  //{vx0,mvx0} for each rp; vx=vx0+mvx*csi
406  double avx[][2]={{0.515483,-1.0123},{0.494122,-1.0534},{0.2217,-1.483},{0.004633,-1.0719}};
407  //{vy0,mvy0} for each rp; vy=vy0+mvy*csi
408  double avy[][2]={{0.371418,-1.6327},{0.349035,-1.6955},{0.0815,-2.59},{0.007592,-4.0841}};
409  //{D0,md,a,b} for each rp; D=D0+(md+a*thetax)*csi+b*thetax
410  double ddx[][4]= {{-0.082336,-0.092513,112.3436,-82.5029},{-0.086927,-0.097670,114.9513,-82.9835},
411  {-0.092117,-0.0915,180.6236,-82.443},{-0.050470,0.058837,208.1106,20.8198}};
412  // {10sigma_x+0.5mm,10sigma_y+0.5mm}
413  double detlim[][2]={{0,0},{0.0028,0.0021},{0,0},{0.0008,0.0013}};
414  //{rmax,dmax}
415  double pipelim[][2]={{0.026,0.026},{0.04,0.04},{0.0226,0.0177},{0.04,0.04}};
416 
417 
418  for(int j=0; j<no_rp ; j++) {
419  //y=Ly*thetay+vy*y0
420  //x=Lx*thetax+vx*x0-csi*D
421  y_par[j]=ThetaY*(leffy[j][0]+leffy[j][1]*csi)+(avy[j][0]+avy[j][1]*csi)*Y_at_0;
422  x_par[j]=ThetaX*(leffx[j][0]+leffx[j][1]*csi)+(avx[j][0]+avx[j][1]*csi)*X_at_0-
423  csi*(ddx[j][0]+(ddx[j][1]+ddx[j][2]*ThetaX)*csi+ddx[j][3]*ThetaX);
424  }
425 
426 
427  //pass TAN@141
428  if (fabs(y_par[0])<pipelim[0][1] && sqrt((y_par[0]*y_par[0])+(x_par[0]*x_par[0]))<pipelim[0][0]) {
429  //pass 149
430  if ((sqrt((y_par[1]*y_par[1])+(x_par[1]*x_par[1]))<pipelim[1][0]) &&
431  (fabs(y_par[1])>detlim[1][1] || x_par[1]>detlim[1][0])) {
432  accettanza = 1;
433  }
434  }
435 
436 
437  //pass TAN@141
438  if (fabs(y_par[0])<pipelim[0][1] && sqrt((y_par[0])*(y_par[0])+(x_par[0])*(x_par[0]))<pipelim[0][0]) {
439  //pass Q5@198
440  if (fabs(y_par[2])<pipelim[2][1] && sqrt((y_par[2]*y_par[2])+(x_par[2]*x_par[2]))<pipelim[2][0]) {
441  //pass 220
442  if ((sqrt((y_par[3]*y_par[3])+(x_par[3]*x_par[3]))<pipelim[3][0]) &&
443  (fabs(y_par[3])>detlim[3][1] || x_par[3]>detlim[3][0])) {
444  accettanza = 1;
445 
446  PosEvo.setX(1000*x_par[3]);
447  PosEvo.setY(1000*y_par[3]);
448  PosEvo.setZ(1000*rp[3]);
449  if(Pos.z()<vz)PosEvo.setZ(-1000*rp[3]);
450  }
451  }
452 
453  }
454 /*
455  LogDebug("ForwardSim") << "\n"
456  << "ACCETTANZA: "<<accettanza << "\n"
457  << "CSI: "<< csi << "\n"
458  << "Theta_X: " << ThetaX << "\n"
459  << "Theta_Y: " << ThetaY << "\n"
460  << "X_at_0: "<< X_at_0 << "\n"
461  << "Y_at_0: "<< Y_at_0 << "\n"
462  << "x_par[3]: "<< x_par[3] << "\n"
463  << "y_par[3]: "<< y_par[3] << "\n"
464  << "pos " << Pos.x() << " " << Pos.y() << " "
465  << Pos.z() << "\n" << "V "<< vx << " " << vy << " "
466  << vz << "\n"
467 */
468 // --------------
469  return PosEvo;
470 }
T sqrt(T t)
Definition: SSEVec.h:46
int j
Definition: DBlmapReader.cc:9
void TotemSD::PrintAll ( )
virtual

Definition at line 171 of file TotemSD.cc.

References LogDebug, and theHC.

171  {
172  LogDebug("ForwardSim") << "TotemSD: Collection " << theHC->GetName();
173  theHC->PrintAllHits();
174 }
#define LogDebug(id)
TotemG4HitCollection * theHC
Definition: TotemSD.h:99
bool TotemSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *   
)
virtual

Implements SensitiveDetector.

Definition at line 101 of file TotemSD.cc.

References CreateNewHit(), CreateNewHitEvo(), edeposit, GetStepInfo(), HitExists(), NULL, ParentId, ParticleType, and unitID.

101  {
102 
103  if (aStep == NULL) {
104  return true;
105  } else {
106  GetStepInfo(aStep);
107  if (HitExists() == false && edeposit>0.) {
108  CreateNewHit();
109  return true;
110  }
111  if (HitExists() == false && (((unitID==1111 || unitID==2222) &&
112  ParentId==0 && ParticleType==2212))) {
113  CreateNewHitEvo();
114  return true;
115  }
116  }
117  return true;
118 }
int ParentId
Definition: TotemSD.h:124
uint32_t unitID
Definition: TotemSD.h:106
#define NULL
Definition: scimark2.h:8
bool HitExists()
Definition: TotemSD.cc:243
float edeposit
Definition: TotemSD.h:112
short ParticleType
Definition: TotemSD.h:119
void CreateNewHitEvo()
Definition: TotemSD.cc:338
void GetStepInfo(G4Step *aStep)
Definition: TotemSD.cc:202
void CreateNewHit()
Definition: TotemSD.cc:285
void TotemSD::ResetForNewPrimary ( )
private

Definition at line 516 of file TotemSD.cc.

References entrancePoint, hitPoint, incidentEnergy, preStepPoint, and SetToLocal().

Referenced by HitExists().

516  {
517 
519  incidentEnergy = preStepPoint->GetKineticEnergy();
520 }
G4ThreeVector entrancePoint
Definition: TotemSD.h:93
float incidentEnergy
Definition: TotemSD.h:94
G4StepPoint * preStepPoint
Definition: TotemSD.h:110
G4ThreeVector SetToLocal(G4ThreeVector globalPoint)
Definition: TotemSD.cc:194
G4ThreeVector hitPoint
Definition: TotemSD.h:113
uint32_t TotemSD::setDetUnitId ( G4Step *  aStep)
virtual

Implements SensitiveDetector.

Definition at line 120 of file TotemSD.cc.

References TotemVDetectorOrganization::GetUnitID(), and numberingScheme.

Referenced by GetStepInfo().

120  {
121 
122  return (numberingScheme == 0 ? 0 : numberingScheme->GetUnitID(aStep));
123 }
virtual uint32_t GetUnitID(const G4Step *aStep) const =0
TotemVDetectorOrganization * numberingScheme
Definition: TotemSD.h:86
G4ThreeVector TotemSD::SetToLocal ( G4ThreeVector  globalPoint)
private

Definition at line 194 of file TotemSD.cc.

References preStepPoint.

Referenced by ResetForNewPrimary().

194  {
195 
196  G4ThreeVector localPoint;
197  const G4VTouchable* touch= preStepPoint->GetTouchable();
198  localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
199  return localPoint;
200 }
G4StepPoint * preStepPoint
Definition: TotemSD.h:110
void TotemSD::StoreHit ( TotemG4Hit hit)
private

Definition at line 505 of file TotemSD.cc.

References primID, and theHC.

Referenced by CreateNewHit(), and CreateNewHitEvo().

505  {
506 
507  if (primID<0) return;
508  if (hit == 0 ) {
509  edm::LogWarning("ForwardSim") << "TotemSD: hit to be stored is NULL !!";
510  return;
511  }
512 
513  theHC->insert( hit );
514 }
TotemG4HitCollection * theHC
Definition: TotemSD.h:99
G4int primID
Definition: TotemSD.h:95
void TotemSD::Summarize ( )
private

Definition at line 522 of file TotemSD.cc.

Referenced by EndOfEvent().

522  {
523 }
void TotemSD::update ( const BeginOfEvent )
privatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent * >.

Definition at line 180 of file TotemSD.cc.

References clearHits(), eventno, and LogDebug.

Referenced by progressbar.ProgressBar::__next__(), relval_steps.Matrix::__setitem__(), relval_steps.Steps::__setitem__(), Vispa.Gui.VispaWidget.VispaWidget::autosize(), Vispa.Views.LineDecayView.LineDecayContainer::createObject(), Vispa.Views.LineDecayView.LineDecayContainer::deselectAllObjects(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::deselectAllWidgets(), Vispa.Gui.VispaWidget.VispaWidget::enableAutosizing(), progressbar.ProgressBar::finish(), Vispa.Gui.MenuWidget.MenuWidget::leaveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseMoveEvent(), Vispa.Gui.MenuWidget.MenuWidget::mouseMoveEvent(), Vispa.Views.LineDecayView.LineDecayContainer::mouseMoveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseReleaseEvent(), Vispa.Views.LineDecayView.LineDecayContainer::objectMoved(), relval_steps.Steps::overwrite(), Vispa.Views.LineDecayView.LineDecayContainer::removeObject(), Vispa.Gui.ConnectableWidget.ConnectableWidget::removePorts(), Vispa.Gui.FindDialog.FindDialog::reset(), Vispa.Gui.PortConnection.PointToPointConnection::select(), Vispa.Gui.VispaWidget.VispaWidget::select(), Vispa.Views.LineDecayView.LineDecayContainer::select(), Vispa.Gui.VispaWidget.VispaWidget::setText(), Vispa.Gui.VispaWidget.VispaWidget::setTitle(), Vispa.Gui.ZoomableWidget.ZoomableWidget::setZoom(), Vispa.Views.LineDecayView.LineDecayContainer::setZoom(), and Vispa.Gui.PortConnection.PointToPointConnection::updateConnection().

180  {
181  LogDebug("ForwardSim") << " Dispatched BeginOfEvent for " << GetName()
182  << " !" ;
183  clearHits();
184  eventno = (*i)()->GetEventID();
185 }
#define LogDebug(id)
virtual void clearHits()
Definition: TotemSD.cc:190
int eventno
Definition: TotemSD.h:127
void TotemSD::update ( const ::EndOfEvent )
private

Definition at line 187 of file TotemSD.cc.

Referenced by progressbar.ProgressBar::__next__(), relval_steps.Matrix::__setitem__(), relval_steps.Steps::__setitem__(), Vispa.Gui.VispaWidget.VispaWidget::autosize(), Vispa.Views.LineDecayView.LineDecayContainer::createObject(), Vispa.Views.LineDecayView.LineDecayContainer::deselectAllObjects(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::deselectAllWidgets(), Vispa.Gui.VispaWidget.VispaWidget::enableAutosizing(), progressbar.ProgressBar::finish(), Vispa.Gui.MenuWidget.MenuWidget::leaveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseMoveEvent(), Vispa.Gui.MenuWidget.MenuWidget::mouseMoveEvent(), Vispa.Views.LineDecayView.LineDecayContainer::mouseMoveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseReleaseEvent(), Vispa.Views.LineDecayView.LineDecayContainer::objectMoved(), relval_steps.Steps::overwrite(), Vispa.Views.LineDecayView.LineDecayContainer::removeObject(), Vispa.Gui.ConnectableWidget.ConnectableWidget::removePorts(), Vispa.Gui.FindDialog.FindDialog::reset(), Vispa.Gui.PortConnection.PointToPointConnection::select(), Vispa.Gui.VispaWidget.VispaWidget::select(), Vispa.Views.LineDecayView.LineDecayContainer::select(), Vispa.Gui.VispaWidget.VispaWidget::setText(), Vispa.Gui.VispaWidget.VispaWidget::setTitle(), Vispa.Gui.ZoomableWidget.ZoomableWidget::setZoom(), Vispa.Views.LineDecayView.LineDecayContainer::setZoom(), and Vispa.Gui.PortConnection.PointToPointConnection::updateConnection().

187  {
188 }
void TotemSD::UpdateHit ( )
private

Definition at line 473 of file TotemSD.cc.

References currentHit, Eloss, LogDebug, postStepPoint, previousUnitID, primaryID, primID, TotemG4Hit::setEnergyLoss(), tsID, tSliceID, and unitID.

Referenced by CreateNewHit(), CreateNewHitEvo(), and HitExists().

473  {
474  //
475  if (Eloss > 0.) {
476  // currentHit->addEnergyDeposit(edepositEM,edepositHAD);
477 
478 #ifdef debug
479  LogDebug("ForwardSim") << "G4TotemT1SD updateHit: add eloss " << Eloss
480  << "\nCurrentHit=" << currentHit
481  << ", PostStepPoint="
482  << postStepPoint->GetPosition();
483 #endif
484 
486  }
487  // if(PostStepPoint->GetPhysicalVolume() != CurrentPV){
488  // currentHit->setExitPoint(SetToLocal(postStepPoint->GetPosition()));
489  // Local3DPoint exit=currentHit->exitPoint();
490 /*
491 #ifdef debug
492  LogDebug("ForwardSim") << "G4TotemT1SD updateHit: exit point "
493  << exit.x() << " " << exit.y() << " " << exit.z();
494 // LogDebug("ForwardSim") << "Energy deposit in Unit " << unitID << " em " << edepositEM/MeV
495 // << " hadronic " << edepositHAD/MeV << " MeV";
496 #endif
497 */
498 
499  // buffer for next steps:
500  tsID = tSliceID;
501  primID = primaryID;
503 }
#define LogDebug(id)
uint32_t unitID
Definition: TotemSD.h:106
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:156
int primaryID
Definition: TotemSD.h:107
uint32_t previousUnitID
Definition: TotemSD.h:106
int tSliceID
Definition: TotemSD.h:107
float Eloss
Definition: TotemSD.h:118
G4StepPoint * postStepPoint
Definition: TotemSD.h:111
TotemG4Hit * currentHit
Definition: TotemSD.h:103
G4int primID
Definition: TotemSD.h:95
int tsID
Definition: TotemSD.h:102

Member Data Documentation

TotemG4Hit* TotemSD::currentHit
private

Definition at line 103 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), HitExists(), and UpdateHit().

G4VPhysicalVolume* TotemSD::currentPV
private

Definition at line 105 of file TotemSD.h.

Referenced by CreateNewHit(), and GetStepInfo().

float TotemSD::edeposit
private

Definition at line 112 of file TotemSD.h.

Referenced by GetStepInfo(), and ProcessHits().

float TotemSD::Eloss
private

Definition at line 118 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), GetStepInfo(), and UpdateHit().

G4ThreeVector TotemSD::entrancePoint
private

Definition at line 93 of file TotemSD.h.

Referenced by ResetForNewPrimary().

int TotemSD::eventno
private

Definition at line 127 of file TotemSD.h.

Referenced by update().

G4int TotemSD::hcID
private

Definition at line 98 of file TotemSD.h.

Referenced by Initialize().

G4ThreeVector TotemSD::hitPoint
private

Definition at line 113 of file TotemSD.h.

Referenced by GetStepInfo(), and ResetForNewPrimary().

float TotemSD::incidentEnergy
private

Definition at line 94 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), and ResetForNewPrimary().

std::string TotemSD::name
private
TotemVDetectorOrganization* TotemSD::numberingScheme
private

Definition at line 86 of file TotemSD.h.

Referenced by setDetUnitId(), TotemSD(), and ~TotemSD().

float TotemSD::Pabs
private

Definition at line 116 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().

int TotemSD::ParentId
private

Definition at line 124 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), GetStepInfo(), and ProcessHits().

short TotemSD::ParticleType
private

Definition at line 119 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), GetStepInfo(), and ProcessHits().

float TotemSD::PhiAtEntry
private

Definition at line 122 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().

G4ThreeVector TotemSD::Posizio
private

Definition at line 115 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().

G4StepPoint* TotemSD::postStepPoint
private

Definition at line 111 of file TotemSD.h.

Referenced by GetStepInfo(), and UpdateHit().

G4StepPoint* TotemSD::preStepPoint
private

Definition at line 110 of file TotemSD.h.

Referenced by GetStepInfo(), ResetForNewPrimary(), and SetToLocal().

uint32_t TotemSD::previousUnitID
private

Definition at line 106 of file TotemSD.h.

Referenced by HitExists(), and UpdateHit().

int TotemSD::primaryID
private

Definition at line 107 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), GetStepInfo(), HitExists(), and UpdateHit().

G4int TotemSD::primID
private

Definition at line 95 of file TotemSD.h.

Referenced by HitExists(), Initialize(), StoreHit(), and UpdateHit().

TrackingSlaveSD* TotemSD::slave
private

Definition at line 85 of file TotemSD.h.

Referenced by clearHits(), EndOfEvent(), fillHits(), TotemSD(), and ~TotemSD().

TotemG4HitCollection* TotemSD::theHC
private

Definition at line 99 of file TotemSD.h.

Referenced by EndOfEvent(), HitExists(), Initialize(), PrintAll(), and StoreHit().

const SimTrackManager* TotemSD::theManager
private

Definition at line 100 of file TotemSD.h.

float TotemSD::ThetaAtEntry
private

Definition at line 121 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().

G4Track* TotemSD::theTrack
private

Definition at line 104 of file TotemSD.h.

Referenced by CreateNewHit(), and GetStepInfo().

float TotemSD::Tof
private

Definition at line 117 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().

int TotemSD::tsID
private

Definition at line 102 of file TotemSD.h.

Referenced by HitExists(), Initialize(), and UpdateHit().

double TotemSD::tSlice
private

Definition at line 108 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().

int TotemSD::tSliceID
private

Definition at line 107 of file TotemSD.h.

Referenced by CreateNewHit(), GetStepInfo(), HitExists(), and UpdateHit().

uint32_t TotemSD::unitID
private

Definition at line 106 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), GetStepInfo(), HitExists(), ProcessHits(), and UpdateHit().

float TotemSD::Vx
private

Definition at line 125 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().

float TotemSD::Vy
private

Definition at line 125 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().

float TotemSD::Vz
private

Definition at line 125 of file TotemSD.h.

Referenced by CreateNewHit(), CreateNewHitEvo(), and GetStepInfo().