CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups 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 SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
 ~TotemSD () override
 
- Public Member Functions inherited from SensitiveTkDetector
 SensitiveTkDetector (const std::string &iname, const SensitiveDetectorCatalog &clg)
 
- 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 SensitiveDetectorCatalog &, 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 46 of file TotemSD.h.

Constructor & Destructor Documentation

TotemSD::TotemSD ( const std::string &  name,
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.

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

Definition at line 84 of file TotemSD.cc.

References numberingScheme, and slave.

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

Member Function Documentation

void TotemSD::clearHits ( )
overridevirtual

Implements SensitiveDetector.

Definition at line 164 of file TotemSD.cc.

References TrackingSlaveSD::Initialize(), and slave.

Referenced by update().

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

Definition at line 249 of file TotemSD.cc.

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

249  {
250 #ifdef EDM_ML_DEBUG
251  edm::LogVerbatim("ForwardSim") << "TotemSD CreateNewHit for PV " << currentPV->GetName()
252  << " PVid = " << currentPV->GetCopyNo() << " Unit " << unitID << "\n primary "
253  << primaryID << " time slice " << tSliceID << " For Track " << theTrack->GetTrackID()
254  << " which is a " << theTrack->GetDefinition()->GetParticleName();
255 
256  if (theTrack->GetTrackID() == 1) {
257  edm::LogVerbatim("ForwardSim") << " of energy " << theTrack->GetTotalEnergy();
258  } else {
259  edm::LogVerbatim("ForwardSim") << " daughter of part. " << theTrack->GetParentID();
260  }
261 
262  if (theTrack->GetCreatorProcess() != nullptr)
263  edm::LogVerbatim("ForwardSim") << theTrack->GetCreatorProcess()->GetProcessName();
264  else
265  edm::LogVerbatim("ForwardSim") << "NO process";
266 #endif
267 
268  currentHit = new TotemG4Hit;
273 
280 
281  currentHit->setEntry(Posizio.x(), Posizio.y(), Posizio.z());
282 
284  currentHit->setVx(Vx);
285  currentHit->setVy(Vy);
286  currentHit->setVz(Vz);
287 
288  updateHit();
289 
291 }
void setVz(float p)
Definition: TotemG4Hit.cc:178
Log< level::Info, true > LogVerbatim
int ParentId
Definition: TotemSD.h:114
void setPabs(float e)
Definition: TotemG4Hit.cc:148
uint32_t unitID
Definition: TotemSD.h:96
float Pabs
Definition: TotemSD.h:106
G4ThreeVector Posizio
Definition: TotemSD.h:105
float Tof
Definition: TotemSD.h:107
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:150
void setEntry(double x, double y, double z)
Definition: TotemG4Hit.h:55
double tSlice
Definition: TotemSD.h:98
int primaryID
Definition: TotemSD.h:97
float incidentEnergy
Definition: TotemSD.h:85
void storeHit(TotemG4Hit *)
Definition: TotemSD.cc:441
int tSliceID
Definition: TotemSD.h:97
void setVy(float p)
Definition: TotemG4Hit.cc:175
void setTimeSlice(double d)
Definition: TotemG4Hit.cc:133
float Eloss
Definition: TotemSD.h:108
float Vx
Definition: TotemSD.h:115
void setVx(float p)
Definition: TotemG4Hit.cc:172
void setParentId(int p)
Definition: TotemG4Hit.cc:169
TotemG4Hit * currentHit
Definition: TotemSD.h:93
G4VPhysicalVolume * currentPV
Definition: TotemSD.h:95
void setIncidentEnergy(double e)
Definition: TotemG4Hit.cc:124
short ParticleType
Definition: TotemSD.h:109
void setThetaAtEntry(float t)
Definition: TotemG4Hit.cc:156
float PhiAtEntry
Definition: TotemSD.h:112
void setUnitID(uint32_t i)
Definition: TotemG4Hit.cc:130
void setTrackID(int i)
Definition: TotemG4Hit.cc:127
void setPhiAtEntry(float f)
Definition: TotemG4Hit.cc:157
void updateHit()
Definition: TotemSD.cc:424
void setParticleType(short i)
Definition: TotemG4Hit.cc:151
void setTof(float e)
Definition: TotemG4Hit.cc:149
float Vz
Definition: TotemSD.h:115
G4Track * theTrack
Definition: TotemSD.h:94
float Vy
Definition: TotemSD.h:115
float ThetaAtEntry
Definition: TotemSD.h:111
void TotemSD::createNewHitEvo ( )
private

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

293  {
294  // edm::LogVerbatim("ForwardSim") << "INSIDE CREATE NEW HIT EVO ";
295 
296  currentHit = new TotemG4Hit;
301 
308 
309  // edm::LogVerbatim("ForwardSim") << Posizio.x() << " " << Posizio.y() << " " << Posizio.z();
310 
312  currentHit->setVx(Vx);
313  currentHit->setVy(Vy);
314  currentHit->setVz(Vz);
315 
316  G4ThreeVector _PosizioEvo;
317  int flagAcc = 0;
318  _PosizioEvo = posizioEvo(Posizio, Vx, Vy, Vz, Pabs, flagAcc);
319 
320  if (flagAcc == 1) {
321  currentHit->setEntry(_PosizioEvo.x(), _PosizioEvo.y(), _PosizioEvo.z());
322 
323  updateHit();
324 
326  }
327  // edm::LogVerbatim("ForwardSim") << "STORED HIT IN: " << unitID;
328 }
void setVz(float p)
Definition: TotemG4Hit.cc:178
int ParentId
Definition: TotemSD.h:114
void setPabs(float e)
Definition: TotemG4Hit.cc:148
uint32_t unitID
Definition: TotemSD.h:96
float Pabs
Definition: TotemSD.h:106
G4ThreeVector posizioEvo(const G4ThreeVector &, double, double, double, double, int &)
Definition: TotemSD.cc:330
G4ThreeVector Posizio
Definition: TotemSD.h:105
float Tof
Definition: TotemSD.h:107
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:150
void setEntry(double x, double y, double z)
Definition: TotemG4Hit.h:55
double tSlice
Definition: TotemSD.h:98
int primaryID
Definition: TotemSD.h:97
float incidentEnergy
Definition: TotemSD.h:85
void storeHit(TotemG4Hit *)
Definition: TotemSD.cc:441
void setVy(float p)
Definition: TotemG4Hit.cc:175
void setTimeSlice(double d)
Definition: TotemG4Hit.cc:133
float Eloss
Definition: TotemSD.h:108
float Vx
Definition: TotemSD.h:115
void setVx(float p)
Definition: TotemG4Hit.cc:172
void setParentId(int p)
Definition: TotemG4Hit.cc:169
TotemG4Hit * currentHit
Definition: TotemSD.h:93
void setIncidentEnergy(double e)
Definition: TotemG4Hit.cc:124
short ParticleType
Definition: TotemSD.h:109
void setThetaAtEntry(float t)
Definition: TotemG4Hit.cc:156
float PhiAtEntry
Definition: TotemSD.h:112
void setUnitID(uint32_t i)
Definition: TotemG4Hit.cc:130
void setTrackID(int i)
Definition: TotemG4Hit.cc:127
void setPhiAtEntry(float f)
Definition: TotemG4Hit.cc:157
void updateHit()
Definition: TotemSD.cc:424
void setParticleType(short i)
Definition: TotemG4Hit.cc:151
void setTof(float e)
Definition: TotemG4Hit.cc:149
float Vz
Definition: TotemSD.h:115
float Vy
Definition: TotemSD.h:115
float ThetaAtEntry
Definition: TotemSD.h:111
void TotemSD::EndOfEvent ( G4HCofThisEvent *  eventHC)
override

Definition at line 119 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, TrackingSlaveSD::processHits(), slave, and theHC.

119  {
120  // here we loop over transient hits and make them persistent
121  int thehc_entries = theHC->entries();
122  for (int j = 0; j < thehc_entries && j < 15000; j++) {
123  TotemG4Hit* aHit = (*theHC)[j];
124 #ifdef EDM_ML_DEBUG
125  edm::LogVerbatim("ForwardSim") << "HIT NUMERO " << j << "unit ID = " << aHit->getUnitID()
126  << "\n enrty z " << aHit->getEntry().z() << "\n theta "
127  << 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 }
Log< level::Info, true > LogVerbatim
TotemG4HitCollection * theHC
Definition: TotemSD.h:89
float getTof() const
Definition: TotemG4Hit.cc:144
float getEnergyLoss() const
Definition: TotemG4Hit.cc:145
math::XYZPoint getEntry() const
Definition: TotemG4Hit.cc:115
float getPabs() const
Definition: TotemG4Hit.cc:143
float getPhiAtEntry() const
Definition: TotemG4Hit.cc:154
int getTrackID() const
Definition: TotemG4Hit.cc:126
uint32_t getUnitID() const
Definition: TotemG4Hit.cc:129
TrackingSlaveSD * slave
Definition: TotemSD.h:76
int getParticleType() const
Definition: TotemG4Hit.cc:146
virtual bool processHits(const PSimHit &)
float getThetaAtEntry() const
Definition: TotemG4Hit.cc:153
void TotemSD::fillHits ( edm::PSimHitContainer cc,
const std::string &  hname 
)
overridevirtual

Implements SensitiveTkDetector.

Definition at line 151 of file TotemSD.cc.

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

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

Definition at line 173 of file TotemSD.cc.

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

Referenced by ProcessHits().

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

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

210  {
211  if (primaryID < 1) {
212  edm::LogWarning("ForwardSim") << "***** TotemSD error: primaryID = " << primaryID << " maybe detector name changed";
213  }
214 
215  // Update if in the same detector, time-slice and for same track
216  // if (primaryID == primID && tSliceID == tsID && unitID==previousUnitID) {
217  if (tSliceID == tsID && unitID == previousUnitID) {
218  updateHit();
219  return true;
220  }
221 
222  // Reset entry point for new primary
223  if (primaryID != primID)
225 
226  //look in the HitContainer whether a hit with the same primID, unitID,
227  //tSliceID already exists:
228 
229  bool found = false;
230  int thehc_entries = theHC->entries();
231  for (int j = 0; j < thehc_entries && !found; j++) {
232  TotemG4Hit* aPreviousHit = (*theHC)[j];
233  if (aPreviousHit->getTrackID() == primaryID && aPreviousHit->getTimeSliceID() == tSliceID &&
234  aPreviousHit->getUnitID() == unitID) {
235  currentHit = aPreviousHit;
236  found = true;
237  break;
238  }
239  }
240 
241  if (found) {
242  updateHit();
243  return true;
244  } else {
245  return false;
246  }
247 }
void resetForNewPrimary()
Definition: TotemSD.cc:452
uint32_t unitID
Definition: TotemSD.h:96
TotemG4HitCollection * theHC
Definition: TotemSD.h:89
int primaryID
Definition: TotemSD.h:97
uint32_t previousUnitID
Definition: TotemSD.h:96
int tSliceID
Definition: TotemSD.h:97
int getTrackID() const
Definition: TotemG4Hit.cc:126
int getTimeSliceID() const
Definition: TotemG4Hit.cc:134
TotemG4Hit * currentHit
Definition: TotemSD.h:93
uint32_t getUnitID() const
Definition: TotemG4Hit.cc:129
void updateHit()
Definition: TotemSD.cc:424
G4int primID
Definition: TotemSD.h:86
Log< level::Warning, false > LogWarning
int tsID
Definition: TotemSD.h:92
void TotemSD::Initialize ( G4HCofThisEvent *  HCE)
override

Definition at line 105 of file TotemSD.cc.

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

105  {
106 #ifdef EDM_ML_DEBUG
107  edm::LogVerbatim("ForwardSim") << "TotemSD : Initialize called for " << GetName();
108 #endif
109 
110  theHC = new TotemG4HitCollection(GetName(), collectionName[0]);
111  if (hcID < 0)
112  hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
113  HCE->AddHitsCollection(hcID, theHC);
114 
115  tsID = -2;
116  primID = -2;
117 }
Log< level::Info, true > LogVerbatim
TotemG4HitCollection * theHC
Definition: TotemSD.h:89
G4int hcID
Definition: TotemSD.h:88
std::string const collectionName[nCollections]
Definition: Collections.h:59
G4int primID
Definition: TotemSD.h:86
int tsID
Definition: TotemSD.h:92
G4THitsCollection< TotemG4Hit > TotemG4HitCollection
G4ThreeVector TotemSD::posizioEvo ( const G4ThreeVector &  Pos,
double  vx,
double  vy,
double  vz,
double  pabs,
int &  accettanza 
)
private

Definition at line 330 of file TotemSD.cc.

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

Referenced by createNewHitEvo().

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

144  {
145 #ifdef EDM_ML_DEBUG
146  edm::LogVerbatim("ForwardSim") << "TotemSD: Collection " << theHC->GetName();
147 #endif
148  theHC->PrintAllHits();
149 }
Log< level::Info, true > LogVerbatim
TotemG4HitCollection * theHC
Definition: TotemSD.h:89
bool TotemSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *   
)
overridevirtual

Implements SensitiveDetector.

Definition at line 89 of file TotemSD.cc.

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

89  {
90  getStepInfo(aStep);
91  if (!hitExists() && edeposit > 0.) {
92  createNewHit();
93  return true;
94  }
95  if (!hitExists() && (((unitID == 1111 || unitID == 2222) && ParentId == 0 && ParticleType == 2212))) {
97  }
98  return true;
99 }
int ParentId
Definition: TotemSD.h:114
uint32_t unitID
Definition: TotemSD.h:96
void createNewHit()
Definition: TotemSD.cc:249
float edeposit
Definition: TotemSD.h:102
short ParticleType
Definition: TotemSD.h:109
void createNewHitEvo()
Definition: TotemSD.cc:293
bool hitExists()
Definition: TotemSD.cc:210
void getStepInfo(const G4Step *aStep)
Definition: TotemSD.cc:173
void TotemSD::resetForNewPrimary ( )
private

Definition at line 452 of file TotemSD.cc.

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

Referenced by hitExists().

452  {
454  incidentEnergy = preStepPoint->GetKineticEnergy();
455 }
G4ThreeVector setToLocal(const G4ThreeVector &globalPoint)
Definition: TotemSD.cc:166
G4ThreeVector entrancePoint
Definition: TotemSD.h:84
float incidentEnergy
Definition: TotemSD.h:85
const G4StepPoint * preStepPoint
Definition: TotemSD.h:100
G4ThreeVector hitPoint
Definition: TotemSD.h:103
uint32_t TotemSD::setDetUnitId ( const G4Step *  aStep)
overridevirtual

Implements SensitiveDetector.

Definition at line 101 of file TotemSD.cc.

References TotemVDetectorOrganization::getUnitID(), and numberingScheme.

Referenced by getStepInfo().

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

Definition at line 166 of file TotemSD.cc.

References preStepPoint.

Referenced by resetForNewPrimary().

166  {
167  G4ThreeVector localPoint;
168  const G4VTouchable* touch = preStepPoint->GetTouchable();
169  localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
170  return localPoint;
171 }
const G4StepPoint * preStepPoint
Definition: TotemSD.h:100
void TotemSD::storeHit ( TotemG4Hit hit)
private

Definition at line 441 of file TotemSD.cc.

References primID, and theHC.

Referenced by createNewHit(), and createNewHitEvo().

441  {
442  if (primID < 0)
443  return;
444  if (hit == nullptr) {
445  edm::LogWarning("ForwardSim") << "TotemSD: hit to be stored is NULL !!";
446  return;
447  }
448 
449  theHC->insert(hit);
450 }
TotemG4HitCollection * theHC
Definition: TotemSD.h:89
G4int primID
Definition: TotemSD.h:86
Log< level::Warning, false > LogWarning
void TotemSD::update ( const BeginOfEvent )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent * >.

Definition at line 157 of file TotemSD.cc.

References clearHits().

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

157  {
158 #ifdef EDM_ML_DEBUG
159  edm::LogVerbatim("ForwardSim") << " Dispatched BeginOfEvent for " << GetName() << " !";
160 #endif
161  clearHits();
162 }
Log< level::Info, true > LogVerbatim
void clearHits() override
Definition: TotemSD.cc:164
void TotemSD::updateHit ( )
private

Definition at line 424 of file TotemSD.cc.

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

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

424  {
425  //
426  if (Eloss > 0.) {
427 #ifdef EDM_ML_DEBUG
428  edm::LogVerbatim("ForwardSim") << "G4TotemT1SD updateHit: add eloss " << Eloss << "\nCurrentHit=" << currentHit
429  << ", PostStepPoint=" << postStepPoint->GetPosition();
430 #endif
431 
433  }
434 
435  // buffer for next steps:
436  tsID = tSliceID;
437  primID = primaryID;
439 }
Log< level::Info, true > LogVerbatim
uint32_t unitID
Definition: TotemSD.h:96
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:150
int primaryID
Definition: TotemSD.h:97
uint32_t previousUnitID
Definition: TotemSD.h:96
int tSliceID
Definition: TotemSD.h:97
float Eloss
Definition: TotemSD.h:108
TotemG4Hit * currentHit
Definition: TotemSD.h:93
const G4StepPoint * postStepPoint
Definition: TotemSD.h:101
G4int primID
Definition: TotemSD.h:86
int tsID
Definition: TotemSD.h:92

Member Data Documentation

TotemG4Hit* TotemSD::currentHit
private

Definition at line 93 of file TotemSD.h.

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

G4VPhysicalVolume* TotemSD::currentPV
private

Definition at line 95 of file TotemSD.h.

Referenced by createNewHit(), and getStepInfo().

float TotemSD::edeposit
private

Definition at line 102 of file TotemSD.h.

Referenced by getStepInfo(), and ProcessHits().

float TotemSD::Eloss
private

Definition at line 108 of file TotemSD.h.

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

G4ThreeVector TotemSD::entrancePoint
private

Definition at line 84 of file TotemSD.h.

Referenced by resetForNewPrimary().

G4int TotemSD::hcID
private

Definition at line 88 of file TotemSD.h.

Referenced by Initialize().

G4ThreeVector TotemSD::hitPoint
private

Definition at line 103 of file TotemSD.h.

Referenced by getStepInfo(), and resetForNewPrimary().

float TotemSD::incidentEnergy
private

Definition at line 85 of file TotemSD.h.

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

TotemVDetectorOrganization* TotemSD::numberingScheme
private

Definition at line 77 of file TotemSD.h.

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

float TotemSD::Pabs
private

Definition at line 106 of file TotemSD.h.

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

int TotemSD::ParentId
private

Definition at line 114 of file TotemSD.h.

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

short TotemSD::ParticleType
private

Definition at line 109 of file TotemSD.h.

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

float TotemSD::PhiAtEntry
private

Definition at line 112 of file TotemSD.h.

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

G4ThreeVector TotemSD::Posizio
private

Definition at line 105 of file TotemSD.h.

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

const G4StepPoint* TotemSD::postStepPoint
private

Definition at line 101 of file TotemSD.h.

Referenced by getStepInfo(), and updateHit().

const G4StepPoint* TotemSD::preStepPoint
private

Definition at line 100 of file TotemSD.h.

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

uint32_t TotemSD::previousUnitID
private

Definition at line 96 of file TotemSD.h.

Referenced by hitExists(), and updateHit().

int TotemSD::primaryID
private

Definition at line 97 of file TotemSD.h.

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

G4int TotemSD::primID
private

Definition at line 86 of file TotemSD.h.

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

TrackingSlaveSD* TotemSD::slave
private

Definition at line 76 of file TotemSD.h.

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

TotemG4HitCollection* TotemSD::theHC
private

Definition at line 89 of file TotemSD.h.

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

const SimTrackManager* TotemSD::theManager
private

Definition at line 90 of file TotemSD.h.

float TotemSD::ThetaAtEntry
private

Definition at line 111 of file TotemSD.h.

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

G4Track* TotemSD::theTrack
private

Definition at line 94 of file TotemSD.h.

Referenced by createNewHit(), and getStepInfo().

float TotemSD::Tof
private

Definition at line 107 of file TotemSD.h.

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

int TotemSD::tsID
private

Definition at line 92 of file TotemSD.h.

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

double TotemSD::tSlice
private

Definition at line 98 of file TotemSD.h.

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

int TotemSD::tSliceID
private

Definition at line 97 of file TotemSD.h.

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

uint32_t TotemSD::unitID
private

Definition at line 96 of file TotemSD.h.

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

float TotemSD::Vx
private

Definition at line 115 of file TotemSD.h.

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

float TotemSD::Vy
private

Definition at line 115 of file TotemSD.h.

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

float TotemSD::Vz
private

Definition at line 115 of file TotemSD.h.

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