CMS 3D CMS Logo

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

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

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

Public Member Functions

void clearHits () override
 
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
void fillHits (edm::PSimHitContainer &, const std::string &) override
 
void Initialize (G4HCofThisEvent *HCE) override
 
void PrintAll () override
 
bool ProcessHits (G4Step *, G4TouchableHistory *) override
 
uint32_t setDetUnitId (const G4Step *) override
 
 TotemSD (const std::string &, const edm::EventSetup &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
 ~TotemSD () override
 
- Public Member Functions inherited from SensitiveTkDetector
 SensitiveTkDetector (const std::string &iname, const edm::EventSetup &es, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
 
- Public Member Functions inherited from SensitiveDetector
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
const std::vector< std::string > & getNames () const
 
void Initialize (G4HCofThisEvent *eventHC) override
 
bool isCaloSD () const
 
 SensitiveDetector (const std::string &iname, const edm::EventSetup &es, const SensitiveDetectorCatalog &, edm::ParameterSet const &p, bool calo)
 
 ~SensitiveDetector () override
 
- Public Member Functions inherited from Observer< const BeginOfEvent * >
 Observer ()
 
void slotForUpdate (const BeginOfEvent * iT)
 
virtual ~Observer ()
 

Protected Member Functions

void update (const BeginOfEvent *) override
 This routine will be called when the appropriate signal arrives. More...
 
- Protected Member Functions inherited from SensitiveDetector
TrackInformationcmsTrackInformation (const G4Track *aTrack)
 
Local3DPoint ConvertToLocal3DPoint (const G4ThreeVector &point) const
 
Local3DPoint FinalStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint InitialStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint LocalPostStepPosition (const G4Step *step) const
 
Local3DPoint LocalPreStepPosition (const G4Step *step) const
 
void NaNTrap (const G4Step *step) const
 
void setNames (const std::vector< std::string > &)
 

Private Member Functions

void createNewHit ()
 
void createNewHitEvo ()
 
void getStepInfo (const G4Step *aStep)
 
bool hitExists ()
 
G4ThreeVector posizioEvo (const G4ThreeVector &, double, double, double, double, int &)
 
void resetForNewPrimary ()
 
G4ThreeVector setToLocal (const G4ThreeVector &globalPoint)
 
void storeHit (TotemG4Hit *)
 
void updateHit ()
 

Private Attributes

TotemG4HitcurrentHit
 
G4VPhysicalVolume * currentPV
 
float edeposit
 
float Eloss
 
G4ThreeVector entrancePoint
 
G4int hcID
 
G4ThreeVector hitPoint
 
float incidentEnergy
 
TotemVDetectorOrganizationnumberingScheme
 
float Pabs
 
int ParentId
 
short ParticleType
 
float PhiAtEntry
 
G4ThreeVector Posizio
 
const G4StepPoint * postStepPoint
 
const 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

- Protected Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 

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 ( const std::string &  name,
const edm::EventSetup es,
const SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)

Definition at line 47 of file TotemSD.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), numberingScheme, and slave.

52  : SensitiveTkDetector(name, es, clg, p),
53  numberingScheme(nullptr),
54  hcID(-1),
55  theHC(nullptr),
56  theManager(manager),
57  currentHit(nullptr),
58  theTrack(nullptr),
59  currentPV(nullptr),
60  unitID(0),
61  previousUnitID(0),
62  preStepPoint(nullptr),
63  postStepPoint(nullptr) {
64  //Parameters
65  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("TotemSD");
66  int verbn = m_p.getUntrackedParameter<int>("Verbosity");
67 
68  SetVerboseLevel(verbn);
69 
70  slave = new TrackingSlaveSD(name);
71 
72  if (name == "TotemHitsT1") {
74  } else if (name == "TotemHitsT2Si") {
76  } else if (name == "TotemHitsT2Gem") {
78  } else if (name == "TotemHitsRP") {
80  } else {
81  edm::LogWarning("ForwardSim") << "TotemSD: ReadoutName not supported\n";
82  }
83 }
T getUntrackedParameter(std::string const &, T const &) const
uint32_t unitID
Definition: TotemSD.h:97
TotemG4HitCollection * theHC
Definition: TotemSD.h:90
SensitiveTkDetector(const std::string &iname, const edm::EventSetup &es, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
uint32_t previousUnitID
Definition: TotemSD.h:97
TotemVDetectorOrganization * numberingScheme
Definition: TotemSD.h:78
G4int hcID
Definition: TotemSD.h:89
const G4StepPoint * preStepPoint
Definition: TotemSD.h:101
TotemG4Hit * currentHit
Definition: TotemSD.h:94
G4VPhysicalVolume * currentPV
Definition: TotemSD.h:96
const G4StepPoint * postStepPoint
Definition: TotemSD.h:102
TrackingSlaveSD * slave
Definition: TotemSD.h:77
G4Track * theTrack
Definition: TotemSD.h:95
const SimTrackManager * theManager
Definition: TotemSD.h:91
TotemSD::~TotemSD ( )
override

Definition at line 85 of file TotemSD.cc.

References numberingScheme, and slave.

85  {
86  delete slave;
87  delete numberingScheme;
88 }
TotemVDetectorOrganization * numberingScheme
Definition: TotemSD.h:78
TrackingSlaveSD * slave
Definition: TotemSD.h:77

Member Function Documentation

void TotemSD::clearHits ( )
overridevirtual

Implements SensitiveDetector.

Definition at line 160 of file TotemSD.cc.

References TrackingSlaveSD::Initialize(), and slave.

Referenced by update().

160 { slave->Initialize(); }
virtual void Initialize()
TrackingSlaveSD * slave
Definition: TotemSD.h:77
void TotemSD::createNewHit ( )
private

Definition at line 245 of file TotemSD.cc.

References currentHit, currentPV, Eloss, incidentEnergy, LogDebug, 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().

245  {
246 #ifdef debug
247  LogDebug("ForwardSim") << "TotemSD CreateNewHit for"
248  << " PV " << currentPV->GetName() << " PVid = " << currentPV->GetCopyNo()
249  << " MVid = " << currentPV->GetMother()->GetCopyNo() << " Unit " << unitID << "\n"
250  << " primary " << primaryID << " time slice " << tSliceID << " For Track "
251  << theTrack->GetTrackID() << " which is a " << theTrack->GetDefinition()->GetParticleName();
252 
253  if (theTrack->GetTrackID() == 1) {
254  LogDebug("ForwardSim") << " of energy " << theTrack->GetTotalEnergy();
255  } else {
256  LogDebug("ForwardSim") << " daughter of part. " << theTrack->GetParentID();
257  }
258 
259  if (theTrack->GetCreatorProcess() != nullptr)
260  LogDebug("ForwardSim") << theTrack->GetCreatorProcess()->GetProcessName();
261  else
262  LogDebug("ForwardSim") << "NO process";
263 #endif
264 
265  currentHit = new TotemG4Hit;
270 
277 
278  currentHit->setEntry(Posizio.x(), Posizio.y(), Posizio.z());
279 
281  currentHit->setVx(Vx);
282  currentHit->setVy(Vy);
283  currentHit->setVz(Vz);
284 
285  updateHit();
286 
288 }
#define LogDebug(id)
void setVz(float p)
Definition: TotemG4Hit.cc:177
int ParentId
Definition: TotemSD.h:115
void setPabs(float e)
Definition: TotemG4Hit.cc:147
uint32_t unitID
Definition: TotemSD.h:97
float Pabs
Definition: TotemSD.h:107
G4ThreeVector Posizio
Definition: TotemSD.h:106
float Tof
Definition: TotemSD.h:108
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:149
void setEntry(double x, double y, double z)
Definition: TotemG4Hit.h:55
double tSlice
Definition: TotemSD.h:99
int primaryID
Definition: TotemSD.h:98
float incidentEnergy
Definition: TotemSD.h:86
void storeHit(TotemG4Hit *)
Definition: TotemSD.cc:438
int tSliceID
Definition: TotemSD.h:98
void setVy(float p)
Definition: TotemG4Hit.cc:174
void setTimeSlice(double d)
Definition: TotemG4Hit.cc:132
float Eloss
Definition: TotemSD.h:109
float Vx
Definition: TotemSD.h:116
void setVx(float p)
Definition: TotemG4Hit.cc:171
void setParentId(int p)
Definition: TotemG4Hit.cc:168
TotemG4Hit * currentHit
Definition: TotemSD.h:94
G4VPhysicalVolume * currentPV
Definition: TotemSD.h:96
void setIncidentEnergy(double e)
Definition: TotemG4Hit.cc:123
short ParticleType
Definition: TotemSD.h:110
void setThetaAtEntry(float t)
Definition: TotemG4Hit.cc:155
float PhiAtEntry
Definition: TotemSD.h:113
void setUnitID(uint32_t i)
Definition: TotemG4Hit.cc:129
void setTrackID(int i)
Definition: TotemG4Hit.cc:126
void setPhiAtEntry(float f)
Definition: TotemG4Hit.cc:156
void updateHit()
Definition: TotemSD.cc:421
void setParticleType(short i)
Definition: TotemG4Hit.cc:150
void setTof(float e)
Definition: TotemG4Hit.cc:148
float Vz
Definition: TotemSD.h:116
G4Track * theTrack
Definition: TotemSD.h:95
float Vy
Definition: TotemSD.h:116
float ThetaAtEntry
Definition: TotemSD.h:112
void TotemSD::createNewHitEvo ( )
private

Definition at line 290 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().

290  {
291  // LogDebug("ForwardSim") << "INSIDE CREATE NEW HIT EVO ";
292 
293  currentHit = new TotemG4Hit;
298 
305 
306  // LogDebug("ForwardSim") << Posizio.x() << " " << Posizio.y() << " " << Posizio.z();
307 
309  currentHit->setVx(Vx);
310  currentHit->setVy(Vy);
311  currentHit->setVz(Vz);
312 
313  G4ThreeVector _PosizioEvo;
314  int flagAcc = 0;
315  _PosizioEvo = posizioEvo(Posizio, Vx, Vy, Vz, Pabs, flagAcc);
316 
317  if (flagAcc == 1) {
318  currentHit->setEntry(_PosizioEvo.x(), _PosizioEvo.y(), _PosizioEvo.z());
319 
320  updateHit();
321 
323  }
324  // LogDebug("ForwardSim") << "STORED HIT IN: " << unitID;
325 }
void setVz(float p)
Definition: TotemG4Hit.cc:177
int ParentId
Definition: TotemSD.h:115
void setPabs(float e)
Definition: TotemG4Hit.cc:147
uint32_t unitID
Definition: TotemSD.h:97
float Pabs
Definition: TotemSD.h:107
G4ThreeVector posizioEvo(const G4ThreeVector &, double, double, double, double, int &)
Definition: TotemSD.cc:327
G4ThreeVector Posizio
Definition: TotemSD.h:106
float Tof
Definition: TotemSD.h:108
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:149
void setEntry(double x, double y, double z)
Definition: TotemG4Hit.h:55
double tSlice
Definition: TotemSD.h:99
int primaryID
Definition: TotemSD.h:98
float incidentEnergy
Definition: TotemSD.h:86
void storeHit(TotemG4Hit *)
Definition: TotemSD.cc:438
void setVy(float p)
Definition: TotemG4Hit.cc:174
void setTimeSlice(double d)
Definition: TotemG4Hit.cc:132
float Eloss
Definition: TotemSD.h:109
float Vx
Definition: TotemSD.h:116
void setVx(float p)
Definition: TotemG4Hit.cc:171
void setParentId(int p)
Definition: TotemG4Hit.cc:168
TotemG4Hit * currentHit
Definition: TotemSD.h:94
void setIncidentEnergy(double e)
Definition: TotemG4Hit.cc:123
short ParticleType
Definition: TotemSD.h:110
void setThetaAtEntry(float t)
Definition: TotemG4Hit.cc:155
float PhiAtEntry
Definition: TotemSD.h:113
void setUnitID(uint32_t i)
Definition: TotemG4Hit.cc:129
void setTrackID(int i)
Definition: TotemG4Hit.cc:126
void setPhiAtEntry(float f)
Definition: TotemG4Hit.cc:156
void updateHit()
Definition: TotemSD.cc:421
void setParticleType(short i)
Definition: TotemG4Hit.cc:150
void setTof(float e)
Definition: TotemG4Hit.cc:148
float Vz
Definition: TotemSD.h:116
float Vy
Definition: TotemSD.h:116
float ThetaAtEntry
Definition: TotemSD.h:112
void TotemSD::EndOfEvent ( G4HCofThisEvent *  eventHC)
override

Definition at line 118 of file TotemSD.cc.

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

118  {
119  // here we loop over transient hits and make them persistent
120  for (int j = 0; j < theHC->entries() && j < 15000; j++) {
121  TotemG4Hit* aHit = (*theHC)[j];
122 #ifdef ddebug
123  LogDebug("ForwardSim") << "HIT NUMERO " << j << "unit ID = " << aHit->getUnitID() << "\n"
124  << " "
125  << "enrty z " << aHit->getEntry().z() << "\n"
126  << " "
127  << "theta " << aHit->getThetaAtEntry() << "\n";
128 #endif
129  Local3DPoint theExitPoint(0, 0, 0);
130  Local3DPoint Entrata(aHit->getEntry().x(), aHit->getEntry().y(), aHit->getEntry().z());
131  slave->processHits(PSimHit(Entrata,
132  theExitPoint,
133  aHit->getPabs(),
134  aHit->getTof(),
135  aHit->getEnergyLoss(),
136  aHit->getParticleType(),
137  aHit->getUnitID(),
138  aHit->getTrackID(),
139  aHit->getThetaAtEntry(),
140  aHit->getPhiAtEntry()));
141  }
142 }
#define LogDebug(id)
TotemG4HitCollection * theHC
Definition: TotemSD.h:90
float getTof() const
Definition: TotemG4Hit.cc:143
float getEnergyLoss() const
Definition: TotemG4Hit.cc:144
math::XYZPoint getEntry() const
Definition: TotemG4Hit.cc:114
float getPabs() const
Definition: TotemG4Hit.cc:142
float getPhiAtEntry() const
Definition: TotemG4Hit.cc:153
int getTrackID() const
Definition: TotemG4Hit.cc:125
uint32_t getUnitID() const
Definition: TotemG4Hit.cc:128
TrackingSlaveSD * slave
Definition: TotemSD.h:77
int getParticleType() const
Definition: TotemG4Hit.cc:145
virtual bool processHits(const PSimHit &)
float getThetaAtEntry() const
Definition: TotemG4Hit.cc:152
void TotemSD::fillHits ( edm::PSimHitContainer cc,
const std::string &  hname 
)
overridevirtual

Implements SensitiveTkDetector.

Definition at line 149 of file TotemSD.cc.

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

149  {
150  if (slave->name() == hname) {
151  cc = slave->hits();
152  }
153 }
std::string name() const
std::vector< PSimHit > & hits()
TrackingSlaveSD * slave
Definition: TotemSD.h:77
void TotemSD::getStepInfo ( const G4Step *  aStep)
private

Definition at line 169 of file TotemSD.cc.

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

Referenced by ProcessHits().

169  {
170  preStepPoint = aStep->GetPreStepPoint();
171  postStepPoint = aStep->GetPostStepPoint();
172  theTrack = aStep->GetTrack();
173  hitPoint = preStepPoint->GetPosition();
174  currentPV = preStepPoint->GetPhysicalVolume();
175 
176  // double weight = 1;
177  G4String name = currentPV->GetName();
178  name.assign(name, 0, 4);
179  G4String particleType = theTrack->GetDefinition()->GetParticleName();
180  edeposit = aStep->GetTotalEnergyDeposit();
181 
182  tSlice = (postStepPoint->GetGlobalTime()) / nanosecond;
183  tSliceID = (int)tSlice;
184  unitID = setDetUnitId(aStep);
185 #ifdef debug
186  LogDebug("ForwardSim") << "UNITa " << unitID;
187 #endif
188  primaryID = theTrack->GetTrackID();
189 
190  Posizio = hitPoint;
191  Pabs = aStep->GetPreStepPoint()->GetMomentum().mag() / GeV;
192  Tof = aStep->GetPostStepPoint()->GetGlobalTime() / nanosecond;
193 
194  Eloss = aStep->GetTotalEnergyDeposit() / GeV;
195  ParticleType = theTrack->GetDefinition()->GetPDGEncoding();
196 
197  ThetaAtEntry = aStep->GetPreStepPoint()->GetPosition().theta() / deg;
198  PhiAtEntry = aStep->GetPreStepPoint()->GetPosition().phi() / deg;
199 
200  ParentId = theTrack->GetParentID();
201  Vx = theTrack->GetVertexPosition().x();
202  Vy = theTrack->GetVertexPosition().y();
203  Vz = theTrack->GetVertexPosition().z();
204 }
#define LogDebug(id)
int ParentId
Definition: TotemSD.h:115
uint32_t unitID
Definition: TotemSD.h:97
const double GeV
Definition: MathUtil.h:16
float Pabs
Definition: TotemSD.h:107
G4ThreeVector Posizio
Definition: TotemSD.h:106
float Tof
Definition: TotemSD.h:108
double tSlice
Definition: TotemSD.h:99
int primaryID
Definition: TotemSD.h:98
int tSliceID
Definition: TotemSD.h:98
const G4StepPoint * preStepPoint
Definition: TotemSD.h:101
float Eloss
Definition: TotemSD.h:109
float edeposit
Definition: TotemSD.h:103
float Vx
Definition: TotemSD.h:116
G4ThreeVector hitPoint
Definition: TotemSD.h:104
G4VPhysicalVolume * currentPV
Definition: TotemSD.h:96
uint32_t setDetUnitId(const G4Step *) override
Definition: TotemSD.cc:102
short ParticleType
Definition: TotemSD.h:110
float PhiAtEntry
Definition: TotemSD.h:113
const G4StepPoint * postStepPoint
Definition: TotemSD.h:102
float Vz
Definition: TotemSD.h:116
G4Track * theTrack
Definition: TotemSD.h:95
float Vy
Definition: TotemSD.h:116
float ThetaAtEntry
Definition: TotemSD.h:112
bool TotemSD::hitExists ( )
private

Definition at line 206 of file TotemSD.cc.

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

Referenced by ProcessHits().

206  {
207  if (primaryID < 1) {
208  edm::LogWarning("ForwardSim") << "***** TotemSD error: primaryID = " << primaryID << " maybe detector name changed";
209  }
210 
211  // Update if in the same detector, time-slice and for same track
212  // if (primaryID == primID && tSliceID == tsID && unitID==previousUnitID) {
213  if (tSliceID == tsID && unitID == previousUnitID) {
214  updateHit();
215  return true;
216  }
217 
218  // Reset entry point for new primary
219  if (primaryID != primID)
221 
222  //look in the HitContainer whether a hit with the same primID, unitID,
223  //tSliceID already exists:
224 
225  bool found = false;
226 
227  for (int j = 0; j < theHC->entries() && !found; j++) {
228  TotemG4Hit* aPreviousHit = (*theHC)[j];
229  if (aPreviousHit->getTrackID() == primaryID && aPreviousHit->getTimeSliceID() == tSliceID &&
230  aPreviousHit->getUnitID() == unitID) {
231  currentHit = aPreviousHit;
232  found = true;
233  break;
234  }
235  }
236 
237  if (found) {
238  updateHit();
239  return true;
240  } else {
241  return false;
242  }
243 }
void resetForNewPrimary()
Definition: TotemSD.cc:449
uint32_t unitID
Definition: TotemSD.h:97
TotemG4HitCollection * theHC
Definition: TotemSD.h:90
int primaryID
Definition: TotemSD.h:98
uint32_t previousUnitID
Definition: TotemSD.h:97
int tSliceID
Definition: TotemSD.h:98
int getTrackID() const
Definition: TotemG4Hit.cc:125
int getTimeSliceID() const
Definition: TotemG4Hit.cc:133
TotemG4Hit * currentHit
Definition: TotemSD.h:94
uint32_t getUnitID() const
Definition: TotemG4Hit.cc:128
void updateHit()
Definition: TotemSD.cc:421
G4int primID
Definition: TotemSD.h:87
int tsID
Definition: TotemSD.h:93
void TotemSD::Initialize ( G4HCofThisEvent *  HCE)
override

Definition at line 106 of file TotemSD.cc.

References bysipixelclustmulteventfilter_cfi::collectionName, hcID, LogDebug, primID, theHC, and tsID.

106  {
107  LogDebug("ForwardSim") << "TotemSD : Initialize called for " << GetName();
108 
109  theHC = new TotemG4HitCollection(GetName(), collectionName[0]);
110  if (hcID < 0)
111  hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
112  HCE->AddHitsCollection(hcID, theHC);
113 
114  tsID = -2;
115  primID = -2;
116 }
#define LogDebug(id)
TotemG4HitCollection * theHC
Definition: TotemSD.h:90
G4int hcID
Definition: TotemSD.h:89
G4int primID
Definition: TotemSD.h:87
int tsID
Definition: TotemSD.h:93
G4THitsCollection< TotemG4Hit > TotemG4HitCollection
G4ThreeVector TotemSD::posizioEvo ( const G4ThreeVector &  Pos,
double  vx,
double  vy,
double  vz,
double  pabs,
int &  accettanza 
)
private

Definition at line 327 of file TotemSD.cc.

References funct::abs(), dqmiolumiharvest::j, and mathSSE::sqrt().

Referenced by createNewHitEvo().

328  {
329  accettanza = 0;
330  //Pos.xyz() in mm
331  G4ThreeVector PosEvo;
332  double ThetaX = std::atan((Pos.x() - vx) / (Pos.z() - vz));
333  double ThetaY = std::atan((Pos.y() - vy) / (Pos.z() - vz));
334  double X_at_0 = (vx - ((Pos.x() - vx) / (Pos.z() - vz)) * vz) / 1000.;
335  double Y_at_0 = (vy - ((Pos.y() - vy) / (Pos.z() - vz)) * vz) / 1000.;
336 
337  // double temp_evo_X;
338  // double temp_evo_Y;
339  // double temp_evo_Z;
340  // temp_evo_Z = std::abs(Pos.z())/Pos.z()*220000.;
341 
342  //csi=-dp/d
343  double csi = std::abs((7000. - pabs) / 7000.);
344 
345  // all in m
346  const int no_rp = 4;
347  double x_par[no_rp + 1];
348  double y_par[no_rp + 1];
349  //rp z position
350  double rp[no_rp] = {141., 149., 198., 220.};
351  //{lx0,mlx} for each rp; Lx=lx0+mlx*csi
352  double leffx[][2] = {{122.5429, -46.9312}, {125.4194, -49.1849}, {152.6, -81.157}, {98.8914, -131.8390}};
353  //{ly0,mly} for each rp; Ly=ly0+mly*csi
354  double leffy[][2] = {{124.2314, -55.4852}, {127.7825, -57.4503}, {179.455, -76.274}, {273.0931, -40.4626}};
355  //{vx0,mvx0} for each rp; vx=vx0+mvx*csi
356  double avx[][2] = {{0.515483, -1.0123}, {0.494122, -1.0534}, {0.2217, -1.483}, {0.004633, -1.0719}};
357  //{vy0,mvy0} for each rp; vy=vy0+mvy*csi
358  double avy[][2] = {{0.371418, -1.6327}, {0.349035, -1.6955}, {0.0815, -2.59}, {0.007592, -4.0841}};
359  //{D0,md,a,b} for each rp; D=D0+(md+a*thetax)*csi+b*thetax
360  double ddx[][4] = {{-0.082336, -0.092513, 112.3436, -82.5029},
361  {-0.086927, -0.097670, 114.9513, -82.9835},
362  {-0.092117, -0.0915, 180.6236, -82.443},
363  {-0.050470, 0.058837, 208.1106, 20.8198}};
364  // {10sigma_x+0.5mm,10sigma_y+0.5mm}
365  double detlim[][2] = {{0, 0}, {0.0028, 0.0021}, {0, 0}, {0.0008, 0.0013}};
366  //{rmax,dmax}
367  double pipelim[][2] = {{0.026, 0.026}, {0.04, 0.04}, {0.0226, 0.0177}, {0.04, 0.04}};
368 
369  for (int j = 0; j < no_rp; j++) {
370  //y=Ly*thetay+vy*y0
371  //x=Lx*thetax+vx*x0-csi*D
372  y_par[j] = ThetaY * (leffy[j][0] + leffy[j][1] * csi) + (avy[j][0] + avy[j][1] * csi) * Y_at_0;
373  x_par[j] = ThetaX * (leffx[j][0] + leffx[j][1] * csi) + (avx[j][0] + avx[j][1] * csi) * X_at_0 -
374  csi * (ddx[j][0] + (ddx[j][1] + ddx[j][2] * ThetaX) * csi + ddx[j][3] * ThetaX);
375  }
376 
377  //pass TAN@141
378  if (std::abs(y_par[0]) < pipelim[0][1] && sqrt((y_par[0] * y_par[0]) + (x_par[0] * x_par[0])) < pipelim[0][0]) {
379  //pass 149
380  if ((sqrt((y_par[1] * y_par[1]) + (x_par[1] * x_par[1])) < pipelim[1][0]) &&
381  (std::abs(y_par[1]) > detlim[1][1] || x_par[1] > detlim[1][0])) {
382  accettanza = 1;
383  }
384  }
385 
386  //pass TAN@141
387  if (std::abs(y_par[0]) < pipelim[0][1] && sqrt((y_par[0]) * (y_par[0]) + (x_par[0]) * (x_par[0])) < pipelim[0][0]) {
388  //pass Q5@198
389  if (std::abs(y_par[2]) < pipelim[2][1] && sqrt((y_par[2] * y_par[2]) + (x_par[2] * x_par[2])) < pipelim[2][0]) {
390  //pass 220
391  if ((sqrt((y_par[3] * y_par[3]) + (x_par[3] * x_par[3])) < pipelim[3][0]) &&
392  (std::abs(y_par[3]) > detlim[3][1] || x_par[3] > detlim[3][0])) {
393  accettanza = 1;
394 
395  PosEvo.setX(1000 * x_par[3]);
396  PosEvo.setY(1000 * y_par[3]);
397  PosEvo.setZ(1000 * rp[3]);
398  if (Pos.z() < vz)
399  PosEvo.setZ(-1000 * rp[3]);
400  }
401  }
402  }
403  /*
404  LogDebug("ForwardSim") << "\n"
405  << "ACCETTANZA: "<<accettanza << "\n"
406  << "CSI: "<< csi << "\n"
407  << "Theta_X: " << ThetaX << "\n"
408  << "Theta_Y: " << ThetaY << "\n"
409  << "X_at_0: "<< X_at_0 << "\n"
410  << "Y_at_0: "<< Y_at_0 << "\n"
411  << "x_par[3]: "<< x_par[3] << "\n"
412  << "y_par[3]: "<< y_par[3] << "\n"
413  << "pos " << Pos.x() << " " << Pos.y() << " "
414  << Pos.z() << "\n" << "V "<< vx << " " << vy << " "
415  << vz << "\n"
416 */
417  // --------------
418  return PosEvo;
419 }
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void TotemSD::PrintAll ( )
override

Definition at line 144 of file TotemSD.cc.

References LogDebug, and theHC.

144  {
145  LogDebug("ForwardSim") << "TotemSD: Collection " << theHC->GetName();
146  theHC->PrintAllHits();
147 }
#define LogDebug(id)
TotemG4HitCollection * theHC
Definition: TotemSD.h:90
bool TotemSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *   
)
overridevirtual

Implements SensitiveDetector.

Definition at line 90 of file TotemSD.cc.

References createNewHit(), createNewHitEvo(), edeposit, getStepInfo(), hitExists(), ParentId, ParticleType, and unitID.

90  {
91  getStepInfo(aStep);
92  if (!hitExists() && edeposit > 0.) {
93  createNewHit();
94  return true;
95  }
96  if (!hitExists() && (((unitID == 1111 || unitID == 2222) && ParentId == 0 && ParticleType == 2212))) {
98  }
99  return true;
100 }
int ParentId
Definition: TotemSD.h:115
uint32_t unitID
Definition: TotemSD.h:97
void createNewHit()
Definition: TotemSD.cc:245
float edeposit
Definition: TotemSD.h:103
short ParticleType
Definition: TotemSD.h:110
void createNewHitEvo()
Definition: TotemSD.cc:290
bool hitExists()
Definition: TotemSD.cc:206
void getStepInfo(const G4Step *aStep)
Definition: TotemSD.cc:169
void TotemSD::resetForNewPrimary ( )
private

Definition at line 449 of file TotemSD.cc.

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

Referenced by hitExists().

449  {
451  incidentEnergy = preStepPoint->GetKineticEnergy();
452 }
G4ThreeVector setToLocal(const G4ThreeVector &globalPoint)
Definition: TotemSD.cc:162
G4ThreeVector entrancePoint
Definition: TotemSD.h:85
float incidentEnergy
Definition: TotemSD.h:86
const G4StepPoint * preStepPoint
Definition: TotemSD.h:101
G4ThreeVector hitPoint
Definition: TotemSD.h:104
uint32_t TotemSD::setDetUnitId ( const G4Step *  aStep)
overridevirtual

Implements SensitiveDetector.

Definition at line 102 of file TotemSD.cc.

References TotemVDetectorOrganization::getUnitID(), and numberingScheme.

Referenced by getStepInfo().

102  {
103  return (numberingScheme == nullptr ? 0 : numberingScheme->getUnitID(aStep));
104 }
TotemVDetectorOrganization * numberingScheme
Definition: TotemSD.h:78
virtual uint32_t getUnitID(const G4Step *aStep) const =0
G4ThreeVector TotemSD::setToLocal ( const G4ThreeVector &  globalPoint)
private

Definition at line 162 of file TotemSD.cc.

References preStepPoint.

Referenced by resetForNewPrimary().

162  {
163  G4ThreeVector localPoint;
164  const G4VTouchable* touch = preStepPoint->GetTouchable();
165  localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
166  return localPoint;
167 }
const G4StepPoint * preStepPoint
Definition: TotemSD.h:101
void TotemSD::storeHit ( TotemG4Hit hit)
private

Definition at line 438 of file TotemSD.cc.

References primID, and theHC.

Referenced by createNewHit(), and createNewHitEvo().

438  {
439  if (primID < 0)
440  return;
441  if (hit == nullptr) {
442  edm::LogWarning("ForwardSim") << "TotemSD: hit to be stored is NULL !!";
443  return;
444  }
445 
446  theHC->insert(hit);
447 }
TotemG4HitCollection * theHC
Definition: TotemSD.h:90
G4int primID
Definition: TotemSD.h:87
void TotemSD::update ( const BeginOfEvent )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent * >.

Definition at line 155 of file TotemSD.cc.

References clearHits(), and LogDebug.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.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(), MatrixUtil.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().

155  {
156  LogDebug("ForwardSim") << " Dispatched BeginOfEvent for " << GetName() << " !";
157  clearHits();
158 }
#define LogDebug(id)
void clearHits() override
Definition: TotemSD.cc:160
void TotemSD::updateHit ( )
private

Definition at line 421 of file TotemSD.cc.

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

Referenced by createNewHit(), createNewHitEvo(), and hitExists().

421  {
422  //
423  if (Eloss > 0.) {
424 #ifdef debug
425  LogDebug("ForwardSim") << "G4TotemT1SD updateHit: add eloss " << Eloss << "\nCurrentHit=" << currentHit
426  << ", PostStepPoint=" << postStepPoint->GetPosition();
427 #endif
428 
430  }
431 
432  // buffer for next steps:
433  tsID = tSliceID;
434  primID = primaryID;
436 }
#define LogDebug(id)
uint32_t unitID
Definition: TotemSD.h:97
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:149
int primaryID
Definition: TotemSD.h:98
uint32_t previousUnitID
Definition: TotemSD.h:97
int tSliceID
Definition: TotemSD.h:98
float Eloss
Definition: TotemSD.h:109
TotemG4Hit * currentHit
Definition: TotemSD.h:94
const G4StepPoint * postStepPoint
Definition: TotemSD.h:102
G4int primID
Definition: TotemSD.h:87
int tsID
Definition: TotemSD.h:93

Member Data Documentation

TotemG4Hit* TotemSD::currentHit
private

Definition at line 94 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), hitExists(), and updateHit().

G4VPhysicalVolume* TotemSD::currentPV
private

Definition at line 96 of file TotemSD.h.

Referenced by createNewHit(), and getStepInfo().

float TotemSD::edeposit
private

Definition at line 103 of file TotemSD.h.

Referenced by getStepInfo(), and ProcessHits().

float TotemSD::Eloss
private

Definition at line 109 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), getStepInfo(), and updateHit().

G4ThreeVector TotemSD::entrancePoint
private

Definition at line 85 of file TotemSD.h.

Referenced by resetForNewPrimary().

G4int TotemSD::hcID
private

Definition at line 89 of file TotemSD.h.

Referenced by Initialize().

G4ThreeVector TotemSD::hitPoint
private

Definition at line 104 of file TotemSD.h.

Referenced by getStepInfo(), and resetForNewPrimary().

float TotemSD::incidentEnergy
private

Definition at line 86 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), and resetForNewPrimary().

TotemVDetectorOrganization* TotemSD::numberingScheme
private

Definition at line 78 of file TotemSD.h.

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

float TotemSD::Pabs
private

Definition at line 107 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), and getStepInfo().

int TotemSD::ParentId
private

Definition at line 115 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), getStepInfo(), and ProcessHits().

short TotemSD::ParticleType
private

Definition at line 110 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), getStepInfo(), and ProcessHits().

float TotemSD::PhiAtEntry
private

Definition at line 113 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), and getStepInfo().

G4ThreeVector TotemSD::Posizio
private

Definition at line 106 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), and getStepInfo().

const G4StepPoint* TotemSD::postStepPoint
private

Definition at line 102 of file TotemSD.h.

Referenced by getStepInfo(), and updateHit().

const G4StepPoint* TotemSD::preStepPoint
private

Definition at line 101 of file TotemSD.h.

Referenced by getStepInfo(), resetForNewPrimary(), and setToLocal().

uint32_t TotemSD::previousUnitID
private

Definition at line 97 of file TotemSD.h.

Referenced by hitExists(), and updateHit().

int TotemSD::primaryID
private

Definition at line 98 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), getStepInfo(), hitExists(), and updateHit().

G4int TotemSD::primID
private

Definition at line 87 of file TotemSD.h.

Referenced by hitExists(), Initialize(), storeHit(), and updateHit().

TrackingSlaveSD* TotemSD::slave
private

Definition at line 77 of file TotemSD.h.

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

TotemG4HitCollection* TotemSD::theHC
private

Definition at line 90 of file TotemSD.h.

Referenced by EndOfEvent(), hitExists(), Initialize(), PrintAll(), and storeHit().

const SimTrackManager* TotemSD::theManager
private

Definition at line 91 of file TotemSD.h.

float TotemSD::ThetaAtEntry
private

Definition at line 112 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), and getStepInfo().

G4Track* TotemSD::theTrack
private

Definition at line 95 of file TotemSD.h.

Referenced by createNewHit(), and getStepInfo().

float TotemSD::Tof
private

Definition at line 108 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), and getStepInfo().

int TotemSD::tsID
private

Definition at line 93 of file TotemSD.h.

Referenced by hitExists(), Initialize(), and updateHit().

double TotemSD::tSlice
private

Definition at line 99 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), and getStepInfo().

int TotemSD::tSliceID
private

Definition at line 98 of file TotemSD.h.

Referenced by createNewHit(), getStepInfo(), hitExists(), and updateHit().

uint32_t TotemSD::unitID
private

Definition at line 97 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), getStepInfo(), hitExists(), ProcessHits(), and updateHit().

float TotemSD::Vx
private

Definition at line 116 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), and getStepInfo().

float TotemSD::Vy
private

Definition at line 116 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), and getStepInfo().

float TotemSD::Vz
private

Definition at line 116 of file TotemSD.h.

Referenced by createNewHit(), createNewHitEvo(), and getStepInfo().