CMS 3D CMS Logo

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 DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
 ~TotemSD () override
 
- Public Member Functions inherited from SensitiveTkDetector
 SensitiveTkDetector (const std::string &iname, const DDCompactView &cpv, 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 DDCompactView &cpv, 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 DDCompactView cpv,
const SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)

Definition at line 48 of file TotemSD.cc.

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

50  :
51  SensitiveTkDetector(name, cpv, clg, p), numberingScheme(nullptr),
52  hcID(-1), theHC(nullptr), theManager(manager), currentHit(nullptr), theTrack(nullptr),
53  currentPV(nullptr), unitID(0), previousUnitID(0), preStepPoint(nullptr),
54  postStepPoint(nullptr){
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 
62  slave = new TrackingSlaveSD(name);
63 
64  if (name == "TotemHitsT1") {
66  } else if (name == "TotemHitsT2Si") {
68  } else if (name == "TotemHitsT2Gem") {
70  } else if (name == "TotemHitsRP") {
72  } else {
73  edm::LogWarning("ForwardSim") << "TotemSD: ReadoutName not supported\n";
74  }
75 }
T getUntrackedParameter(std::string const &, T const &) const
uint32_t unitID
Definition: TotemSD.h:100
TotemG4HitCollection * theHC
Definition: TotemSD.h:93
uint32_t previousUnitID
Definition: TotemSD.h:100
TotemVDetectorOrganization * numberingScheme
Definition: TotemSD.h:81
G4int hcID
Definition: TotemSD.h:92
const G4StepPoint * preStepPoint
Definition: TotemSD.h:104
TotemG4Hit * currentHit
Definition: TotemSD.h:97
G4VPhysicalVolume * currentPV
Definition: TotemSD.h:99
const G4StepPoint * postStepPoint
Definition: TotemSD.h:105
TrackingSlaveSD * slave
Definition: TotemSD.h:80
G4Track * theTrack
Definition: TotemSD.h:98
const SimTrackManager * theManager
Definition: TotemSD.h:94
SensitiveTkDetector(const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
TotemSD::~TotemSD ( )
override

Definition at line 77 of file TotemSD.cc.

References numberingScheme, and slave.

77  {
78  delete slave;
79  delete numberingScheme;
80 }
TotemVDetectorOrganization * numberingScheme
Definition: TotemSD.h:81
TrackingSlaveSD * slave
Definition: TotemSD.h:80

Member Function Documentation

void TotemSD::clearHits ( )
overridevirtual

Implements SensitiveDetector.

Definition at line 155 of file TotemSD.cc.

References TrackingSlaveSD::Initialize(), and slave.

Referenced by update().

155  {
156  slave->Initialize();
157 }
virtual void Initialize()
TrackingSlaveSD * slave
Definition: TotemSD.h:80
void TotemSD::createNewHit ( )
private

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

248  {
249 
250 #ifdef debug
251  LogDebug("ForwardSim") << "TotemSD CreateNewHit for"
252  << " PV " << currentPV->GetName()
253  << " PVid = " << currentPV->GetCopyNo()
254  << " MVid = " << currentPV->GetMother()->GetCopyNo()
255  << " Unit " << unitID << "\n"
256  << " primary " << primaryID
257  << " time slice " << tSliceID
258  << " For Track " << theTrack->GetTrackID()
259  << " which is a "
260  << theTrack->GetDefinition()->GetParticleName();
261 
262  if (theTrack->GetTrackID()==1) {
263  LogDebug("ForwardSim") << " of energy " << theTrack->GetTotalEnergy();
264  } else {
265  LogDebug("ForwardSim") << " daughter of part. " << theTrack->GetParentID();
266  }
267 
268  if (theTrack->GetCreatorProcess()!=nullptr)
269  LogDebug("ForwardSim") << theTrack->GetCreatorProcess()->GetProcessName() ;
270  else
271  LogDebug("ForwardSim") << "NO process";
272 #endif
273 
274 
275  currentHit = new TotemG4Hit;
280 
287 
288  currentHit->setEntry(Posizio.x(),Posizio.y(),Posizio.z());
289 
291  currentHit->setVx(Vx);
292  currentHit->setVy(Vy);
293  currentHit->setVz(Vz);
294 
295  updateHit();
296 
298 }
#define LogDebug(id)
void setVz(float p)
Definition: TotemG4Hit.cc:183
int ParentId
Definition: TotemSD.h:118
void setPabs(float e)
Definition: TotemG4Hit.cc:153
uint32_t unitID
Definition: TotemSD.h:100
float Pabs
Definition: TotemSD.h:110
G4ThreeVector Posizio
Definition: TotemSD.h:109
float Tof
Definition: TotemSD.h:111
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:155
void setEntry(double x, double y, double z)
Definition: TotemG4Hit.h:57
double tSlice
Definition: TotemSD.h:102
int primaryID
Definition: TotemSD.h:101
float incidentEnergy
Definition: TotemSD.h:89
void storeHit(TotemG4Hit *)
Definition: TotemSD.cc:453
int tSliceID
Definition: TotemSD.h:101
void setVy(float p)
Definition: TotemG4Hit.cc:180
void setTimeSlice(double d)
Definition: TotemG4Hit.cc:141
float Eloss
Definition: TotemSD.h:112
float Vx
Definition: TotemSD.h:119
void setVx(float p)
Definition: TotemG4Hit.cc:177
void setParentId(int p)
Definition: TotemG4Hit.cc:174
TotemG4Hit * currentHit
Definition: TotemSD.h:97
G4VPhysicalVolume * currentPV
Definition: TotemSD.h:99
void setIncidentEnergy(double e)
Definition: TotemG4Hit.cc:132
short ParticleType
Definition: TotemSD.h:113
void setThetaAtEntry(float t)
Definition: TotemG4Hit.cc:161
float PhiAtEntry
Definition: TotemSD.h:116
void setUnitID(uint32_t i)
Definition: TotemG4Hit.cc:138
void setTrackID(int i)
Definition: TotemG4Hit.cc:135
void setPhiAtEntry(float f)
Definition: TotemG4Hit.cc:162
void updateHit()
Definition: TotemSD.cc:433
void setParticleType(short i)
Definition: TotemG4Hit.cc:156
void setTof(float e)
Definition: TotemG4Hit.cc:154
float Vz
Definition: TotemSD.h:119
G4Track * theTrack
Definition: TotemSD.h:98
float Vy
Definition: TotemSD.h:119
float ThetaAtEntry
Definition: TotemSD.h:115
void TotemSD::createNewHitEvo ( )
private

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

300  {
301 
302 // LogDebug("ForwardSim") << "INSIDE CREATE NEW HIT EVO ";
303 
304  currentHit = new TotemG4Hit;
309 
316 
317  // LogDebug("ForwardSim") << Posizio.x() << " " << Posizio.y() << " " << Posizio.z();
318 
320  currentHit->setVx(Vx);
321  currentHit->setVy(Vy);
322  currentHit->setVz(Vz);
323 
324  G4ThreeVector _PosizioEvo;
325  int flagAcc=0;
326  _PosizioEvo=posizioEvo(Posizio,Vx,Vy,Vz,Pabs,flagAcc);
327 
328  if(flagAcc==1){
329  currentHit->setEntry(_PosizioEvo.x(),_PosizioEvo.y(),_PosizioEvo.z());
330 
331  updateHit();
332 
334  }
335  // LogDebug("ForwardSim") << "STORED HIT IN: " << unitID;
336 }
void setVz(float p)
Definition: TotemG4Hit.cc:183
int ParentId
Definition: TotemSD.h:118
void setPabs(float e)
Definition: TotemG4Hit.cc:153
uint32_t unitID
Definition: TotemSD.h:100
float Pabs
Definition: TotemSD.h:110
G4ThreeVector posizioEvo(const G4ThreeVector &, double, double, double, double, int &)
Definition: TotemSD.cc:338
G4ThreeVector Posizio
Definition: TotemSD.h:109
float Tof
Definition: TotemSD.h:111
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:155
void setEntry(double x, double y, double z)
Definition: TotemG4Hit.h:57
double tSlice
Definition: TotemSD.h:102
int primaryID
Definition: TotemSD.h:101
float incidentEnergy
Definition: TotemSD.h:89
void storeHit(TotemG4Hit *)
Definition: TotemSD.cc:453
void setVy(float p)
Definition: TotemG4Hit.cc:180
void setTimeSlice(double d)
Definition: TotemG4Hit.cc:141
float Eloss
Definition: TotemSD.h:112
float Vx
Definition: TotemSD.h:119
void setVx(float p)
Definition: TotemG4Hit.cc:177
void setParentId(int p)
Definition: TotemG4Hit.cc:174
TotemG4Hit * currentHit
Definition: TotemSD.h:97
void setIncidentEnergy(double e)
Definition: TotemG4Hit.cc:132
short ParticleType
Definition: TotemSD.h:113
void setThetaAtEntry(float t)
Definition: TotemG4Hit.cc:161
float PhiAtEntry
Definition: TotemSD.h:116
void setUnitID(uint32_t i)
Definition: TotemG4Hit.cc:138
void setTrackID(int i)
Definition: TotemG4Hit.cc:135
void setPhiAtEntry(float f)
Definition: TotemG4Hit.cc:162
void updateHit()
Definition: TotemSD.cc:433
void setParticleType(short i)
Definition: TotemG4Hit.cc:156
void setTof(float e)
Definition: TotemG4Hit.cc:154
float Vz
Definition: TotemSD.h:119
float Vy
Definition: TotemSD.h:119
float ThetaAtEntry
Definition: TotemSD.h:115
void TotemSD::EndOfEvent ( G4HCofThisEvent *  eventHC)
override

Definition at line 114 of file TotemSD.cc.

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

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

Implements SensitiveTkDetector.

Definition at line 145 of file TotemSD.cc.

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

145  {
146  if (slave->name() == hname) { cc=slave->hits(); }
147 }
std::string name() const
std::vector< PSimHit > & hits()
TrackingSlaveSD * slave
Definition: TotemSD.h:80
void TotemSD::getStepInfo ( const G4Step *  aStep)
private

Definition at line 167 of file TotemSD.cc.

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

Referenced by ProcessHits().

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

Definition at line 205 of file TotemSD.cc.

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

Referenced by ProcessHits().

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

Definition at line 101 of file TotemSD.cc.

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

101  {
102 
103  LogDebug("ForwardSim") << "TotemSD : Initialize called for " << GetName();
104 
105  theHC = new TotemG4HitCollection(GetName(), collectionName[0]);
106  if (hcID<0)
107  hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
108  HCE->AddHitsCollection(hcID, theHC);
109 
110  tsID = -2;
111  primID = -2;
112 }
#define LogDebug(id)
TotemG4HitCollection * theHC
Definition: TotemSD.h:93
G4int hcID
Definition: TotemSD.h:92
std::string const collectionName[nCollections]
Definition: Collections.h:47
G4int primID
Definition: TotemSD.h:90
int tsID
Definition: TotemSD.h:96
G4THitsCollection< TotemG4Hit > TotemG4HitCollection
G4ThreeVector TotemSD::posizioEvo ( const G4ThreeVector &  Pos,
double  vx,
double  vy,
double  vz,
double  pabs,
int &  accettanza 
)
private

Definition at line 338 of file TotemSD.cc.

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

Referenced by createNewHitEvo().

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

Definition at line 140 of file TotemSD.cc.

References LogDebug, and theHC.

140  {
141  LogDebug("ForwardSim") << "TotemSD: Collection " << theHC->GetName();
142  theHC->PrintAllHits();
143 }
#define LogDebug(id)
TotemG4HitCollection * theHC
Definition: TotemSD.h:93
bool TotemSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *   
)
overridevirtual

Implements SensitiveDetector.

Definition at line 82 of file TotemSD.cc.

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

82  {
83 
84  getStepInfo(aStep);
85  if (!hitExists() && edeposit>0.) {
86  createNewHit();
87  return true;
88  }
89  if (!hitExists() && (((unitID==1111 || unitID==2222) &&
90  ParentId==0 && ParticleType==2212))) {
92  }
93  return true;
94 }
int ParentId
Definition: TotemSD.h:118
uint32_t unitID
Definition: TotemSD.h:100
void createNewHit()
Definition: TotemSD.cc:248
float edeposit
Definition: TotemSD.h:106
short ParticleType
Definition: TotemSD.h:113
void createNewHitEvo()
Definition: TotemSD.cc:300
bool hitExists()
Definition: TotemSD.cc:205
void getStepInfo(const G4Step *aStep)
Definition: TotemSD.cc:167
void TotemSD::resetForNewPrimary ( )
private

Definition at line 464 of file TotemSD.cc.

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

Referenced by hitExists().

464  {
465 
467  incidentEnergy = preStepPoint->GetKineticEnergy();
468 }
G4ThreeVector setToLocal(const G4ThreeVector &globalPoint)
Definition: TotemSD.cc:159
G4ThreeVector entrancePoint
Definition: TotemSD.h:88
float incidentEnergy
Definition: TotemSD.h:89
const G4StepPoint * preStepPoint
Definition: TotemSD.h:104
G4ThreeVector hitPoint
Definition: TotemSD.h:107
uint32_t TotemSD::setDetUnitId ( const G4Step *  aStep)
overridevirtual

Implements SensitiveDetector.

Definition at line 96 of file TotemSD.cc.

References TotemVDetectorOrganization::getUnitID(), and numberingScheme.

Referenced by getStepInfo().

96  {
97 
98  return (numberingScheme == nullptr ? 0 : numberingScheme->getUnitID(aStep));
99 }
TotemVDetectorOrganization * numberingScheme
Definition: TotemSD.h:81
virtual uint32_t getUnitID(const G4Step *aStep) const =0
G4ThreeVector TotemSD::setToLocal ( const G4ThreeVector &  globalPoint)
private

Definition at line 159 of file TotemSD.cc.

References preStepPoint.

Referenced by resetForNewPrimary().

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

Definition at line 453 of file TotemSD.cc.

References primID, and theHC.

Referenced by createNewHit(), and createNewHitEvo().

453  {
454 
455  if (primID<0) return;
456  if (hit == nullptr ) {
457  edm::LogWarning("ForwardSim") << "TotemSD: hit to be stored is NULL !!";
458  return;
459  }
460 
461  theHC->insert( hit );
462 }
TotemG4HitCollection * theHC
Definition: TotemSD.h:93
G4int primID
Definition: TotemSD.h:90
void TotemSD::update ( const BeginOfEvent )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent * >.

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

149  {
150  LogDebug("ForwardSim") << " Dispatched BeginOfEvent for " << GetName()
151  << " !" ;
152  clearHits();
153 }
#define LogDebug(id)
void clearHits() override
Definition: TotemSD.cc:155
void TotemSD::updateHit ( )
private

Definition at line 433 of file TotemSD.cc.

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

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

433  {
434  //
435  if (Eloss > 0.) {
436 
437 #ifdef debug
438  LogDebug("ForwardSim") << "G4TotemT1SD updateHit: add eloss " << Eloss
439  << "\nCurrentHit=" << currentHit
440  << ", PostStepPoint="
441  << postStepPoint->GetPosition();
442 #endif
443 
445  }
446 
447  // buffer for next steps:
448  tsID = tSliceID;
449  primID = primaryID;
451 }
#define LogDebug(id)
uint32_t unitID
Definition: TotemSD.h:100
void setEnergyLoss(float e)
Definition: TotemG4Hit.cc:155
int primaryID
Definition: TotemSD.h:101
uint32_t previousUnitID
Definition: TotemSD.h:100
int tSliceID
Definition: TotemSD.h:101
float Eloss
Definition: TotemSD.h:112
TotemG4Hit * currentHit
Definition: TotemSD.h:97
const G4StepPoint * postStepPoint
Definition: TotemSD.h:105
G4int primID
Definition: TotemSD.h:90
int tsID
Definition: TotemSD.h:96

Member Data Documentation

TotemG4Hit* TotemSD::currentHit
private

Definition at line 97 of file TotemSD.h.

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

G4VPhysicalVolume* TotemSD::currentPV
private

Definition at line 99 of file TotemSD.h.

Referenced by createNewHit(), and getStepInfo().

float TotemSD::edeposit
private

Definition at line 106 of file TotemSD.h.

Referenced by getStepInfo(), and ProcessHits().

float TotemSD::Eloss
private

Definition at line 112 of file TotemSD.h.

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

G4ThreeVector TotemSD::entrancePoint
private

Definition at line 88 of file TotemSD.h.

Referenced by resetForNewPrimary().

G4int TotemSD::hcID
private

Definition at line 92 of file TotemSD.h.

Referenced by Initialize().

G4ThreeVector TotemSD::hitPoint
private

Definition at line 107 of file TotemSD.h.

Referenced by getStepInfo(), and resetForNewPrimary().

float TotemSD::incidentEnergy
private

Definition at line 89 of file TotemSD.h.

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

TotemVDetectorOrganization* TotemSD::numberingScheme
private

Definition at line 81 of file TotemSD.h.

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

float TotemSD::Pabs
private

Definition at line 110 of file TotemSD.h.

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

int TotemSD::ParentId
private

Definition at line 118 of file TotemSD.h.

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

short TotemSD::ParticleType
private

Definition at line 113 of file TotemSD.h.

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

float TotemSD::PhiAtEntry
private

Definition at line 116 of file TotemSD.h.

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

G4ThreeVector TotemSD::Posizio
private

Definition at line 109 of file TotemSD.h.

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

const G4StepPoint* TotemSD::postStepPoint
private

Definition at line 105 of file TotemSD.h.

Referenced by getStepInfo(), and updateHit().

const G4StepPoint* TotemSD::preStepPoint
private

Definition at line 104 of file TotemSD.h.

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

uint32_t TotemSD::previousUnitID
private

Definition at line 100 of file TotemSD.h.

Referenced by hitExists(), and updateHit().

int TotemSD::primaryID
private

Definition at line 101 of file TotemSD.h.

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

G4int TotemSD::primID
private

Definition at line 90 of file TotemSD.h.

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

TrackingSlaveSD* TotemSD::slave
private

Definition at line 80 of file TotemSD.h.

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

TotemG4HitCollection* TotemSD::theHC
private

Definition at line 93 of file TotemSD.h.

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

const SimTrackManager* TotemSD::theManager
private

Definition at line 94 of file TotemSD.h.

float TotemSD::ThetaAtEntry
private

Definition at line 115 of file TotemSD.h.

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

G4Track* TotemSD::theTrack
private

Definition at line 98 of file TotemSD.h.

Referenced by createNewHit(), and getStepInfo().

float TotemSD::Tof
private

Definition at line 111 of file TotemSD.h.

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

int TotemSD::tsID
private

Definition at line 96 of file TotemSD.h.

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

double TotemSD::tSlice
private

Definition at line 102 of file TotemSD.h.

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

int TotemSD::tSliceID
private

Definition at line 101 of file TotemSD.h.

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

uint32_t TotemSD::unitID
private

Definition at line 100 of file TotemSD.h.

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

float TotemSD::Vx
private

Definition at line 119 of file TotemSD.h.

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

float TotemSD::Vy
private

Definition at line 119 of file TotemSD.h.

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

float TotemSD::Vz
private

Definition at line 119 of file TotemSD.h.

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