CMS 3D CMS Logo

TimingSD.cc
Go to the documentation of this file.
1 // File: TimingSD.cc
3 // Date: 02.2006
4 // Description: Sensitive Detector class for Timing
5 // Modifications:
7 
9 
13 
16 
17 #include "G4Step.hh"
18 #include "G4StepPoint.hh"
19 #include "G4Track.hh"
20 #include "G4VPhysicalVolume.hh"
21 #include "G4SDManager.hh"
22 #include "G4VProcess.hh"
23 
24 #include "G4PhysicalConstants.hh"
25 #include "G4SystemOfUnits.hh"
26 
27 #include <vector>
28 #include <iostream>
29 
30 //#define EDM_ML_DEBUG
31 
32 static const float invgev = 1.0 / CLHEP::GeV;
33 static const double invns = 1.0 / CLHEP::nanosecond;
34 static const double invdeg = 1.0 / CLHEP::deg;
35 
36 //-------------------------------------------------------------------
38  : SensitiveTkDetector(name, clg),
39  theManager(manager),
40  theHC(nullptr),
41  currentHit(nullptr),
42  theTrack(nullptr),
43  preStepPoint(nullptr),
44  postStepPoint(nullptr),
45  unitID(0),
46  previousUnitID(0),
47  primID(-2),
48  hcID(-1),
49  tsID(-2),
50  primaryID(0),
51  tSliceID(-1),
52  timeFactor(1.0),
53  energyCut(1.e+9),
54  energyHistoryCut(1.e+9) {
55  slave = new TrackingSlaveSD(name);
57 }
58 
60  delete slave;
61  delete theEnumerator;
62 }
63 
64 void TimingSD::Initialize(G4HCofThisEvent* HCE) {
65  edm::LogVerbatim("TimingSim") << "TimingSD : Initialize called for " << GetName() << " time slice factor "
66  << timeFactor << "\n MC truth cuts in are " << energyCut / CLHEP::GeV << " GeV and "
67  << energyHistoryCut / CLHEP::GeV << " GeV";
68 
69  theHC = new BscG4HitCollection(GetName(), collectionName[0]);
70  if (hcID < 0) {
71  hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
72  }
73  HCE->AddHitsCollection(hcID, theHC);
74 
75  tsID = -2;
76  primID = -2;
77 }
78 
80  if (val <= 0.0) {
81  return;
82  }
83  timeFactor = val;
84 #ifdef EDM_ML_DEBUG
85  edm::LogVerbatim("TimingSim") << "TimingSD : for " << GetName() << " time slice factor is set to " << timeFactor;
86 #endif
87 }
88 
89 void TimingSD::setCuts(double eCut, double historyCut) {
90  if (eCut > 0.) {
91  energyCut = eCut;
92  }
93  if (historyCut > 0.) {
94  energyHistoryCut = historyCut;
95  }
96 #ifdef EDM_ML_DEBUG
97  edm::LogVerbatim("TimingSim") << "TimingSD : for " << GetName() << " MC truth cuts in are " << energyCut / CLHEP::GeV
98  << " GeV and " << energyHistoryCut / CLHEP::GeV << " GeV";
99 #endif
100 }
101 
102 bool TimingSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
103  edeposit = aStep->GetTotalEnergyDeposit();
104  if (edeposit > 0.f) {
105  getStepInfo(aStep);
106  if (!hitExists(aStep)) {
107  createNewHit(aStep);
108  }
109  }
110  return true;
111 }
112 
113 void TimingSD::getStepInfo(const G4Step* aStep) {
114  preStepPoint = aStep->GetPreStepPoint();
115  postStepPoint = aStep->GetPostStepPoint();
116  hitPointExit = postStepPoint->GetPosition();
118  const G4Track* newTrack = aStep->GetTrack();
119 
120  // neutral particles deliver energy post step
121  // charged particle start deliver energy pre step
122  if (newTrack->GetDefinition()->GetPDGCharge() == 0.0) {
125  tof = (float)(postStepPoint->GetGlobalTime() * invns);
126  } else {
127  hitPoint = preStepPoint->GetPosition();
129  tof = (float)(preStepPoint->GetGlobalTime() * invns);
130  }
131 
132 #ifdef EDM_ML_DEBUG
133  double distGlobal =
134  std::sqrt(std::pow(hitPoint.x() - hitPointExit.x(), 2) + std::pow(hitPoint.y() - hitPointExit.y(), 2) +
135  std::pow(hitPoint.z() - hitPointExit.z(), 2));
136  double distLocal = std::sqrt(std::pow(hitPointLocal.x() - hitPointLocalExit.x(), 2) +
139  LogDebug("TimingSim") << "TimingSD:"
140  << "\n Global entry point: " << hitPoint << "\n Global exit point: " << hitPointExit
141  << "\n Global step length: " << distGlobal << "\n Local entry point: " << hitPointLocal
142  << "\n Local exit point: " << hitPointLocalExit << "\n Local step length: " << distLocal;
143  if (std::fabs(distGlobal - distLocal) > 1.e-6) {
144  LogDebug("TimingSim") << "DIFFERENCE IN DISTANCE \n";
145  }
146 #endif
147 
148  incidentEnergy = preStepPoint->GetKineticEnergy();
149 
150  // should MC truth be saved
151  if (newTrack != theTrack) {
152  theTrack = newTrack;
153  TrackInformation* info = nullptr;
154  if (incidentEnergy > energyCut) {
156  info->storeTrack(true);
157  }
159  if (!info) {
161  }
162  info->putInHistory();
163  }
164  }
165 
166  edeposit *= invgev;
169  edepositHAD = 0.f;
170  } else {
171  edepositEM = 0.f;
173  }
174  // time slice is defined for the entry point
175  tSlice = timeFactor * preStepPoint->GetGlobalTime() * invns;
176  tSliceID = (int)tSlice;
177 
178  unitID = setDetUnitId(aStep);
179  primaryID = theTrack->GetTrackID();
180 }
181 
182 bool TimingSD::hitExists(const G4Step* aStep) {
183  if (!currentHit) {
184  return false;
185  }
186 
187  // Update if in the same detector and time-slice
188  if (tSliceID == tsID && unitID == previousUnitID) {
189  updateHit();
190  return true;
191  }
192 
193  //look in the HitContainer whether a hit with the same primID, unitID,
194  //tSliceID already exists:
195 
196  bool found = false;
197  int thehc_entries = theHC->entries();
198  for (int j = 0; j < thehc_entries; ++j) {
199  BscG4Hit* aHit = (*theHC)[j];
200  if (aHit->getTimeSliceID() == tSliceID && aHit->getUnitID() == unitID) {
201  if (checkHit(aStep, aHit)) {
202  currentHit = aHit;
203  found = true;
204  break;
205  }
206  }
207  }
208  if (found) {
209  updateHit();
210  }
211  return found;
212 }
213 
214 bool TimingSD::checkHit(const G4Step*, BscG4Hit* hit) {
215  // change hit info to fastest primary particle
216  if (tof < hit->getTof()) {
217  hit->setTrackID(primaryID);
218  hit->setIncidentEnergy((float)incidentEnergy);
219  hit->setPabs(float(preStepPoint->GetMomentum().mag()) * invgev);
220  hit->setTof(tof);
221  hit->setParticleType(theTrack->GetDefinition()->GetPDGEncoding());
222 
223  float ThetaAtEntry = (float)(hitPointLocal.theta() * invdeg);
224  float PhiAtEntry = (float)(hitPointLocal.phi() * invdeg);
225 
226  hit->setThetaAtEntry(ThetaAtEntry);
227  hit->setPhiAtEntry(PhiAtEntry);
228 
229  hit->setEntry(hitPoint);
230  hit->setEntryLocalP(hitPointLocal);
231  hit->setExitLocalP(hitPointLocalExit);
232 
233  hit->setParentId(theTrack->GetParentID());
234  hit->setProcessId(theEnumerator->processId(theTrack->GetCreatorProcess()));
235 
236  hit->setVertexPosition(theTrack->GetVertexPosition());
237  }
238  return true;
239 }
240 
242  if (primID < 0)
243  return;
244  if (hit == nullptr) {
245  edm::LogWarning("BscSim") << "BscSD: hit to be stored is NULL !!";
246  return;
247  }
248 
249  theHC->insert(hit);
250 }
251 
252 void TimingSD::createNewHit(const G4Step* aStep) {
253 #ifdef EDM_ML_DEBUG
254  const G4VPhysicalVolume* currentPV = preStepPoint->GetPhysicalVolume();
255  edm::LogVerbatim("TimingSim") << "TimingSD CreateNewHit for " << GetName() << " PV " << currentPV->GetName()
256  << " PVid = " << currentPV->GetCopyNo() << " Unit " << unitID << "\n primary "
257  << primaryID << " Tof(ns)= " << tof << " time slice " << tSliceID
258  << " E(GeV)= " << incidentEnergy << " trackID " << theTrack->GetTrackID() << " "
259  << theTrack->GetDefinition()->GetParticleName() << " parentID "
260  << theTrack->GetParentID();
261 
262  if (theTrack->GetCreatorProcess() != nullptr) {
263  edm::LogVerbatim("TimingSim") << theTrack->GetCreatorProcess()->GetProcessName();
264  } else {
265  edm::LogVerbatim("TimingSim") << " is primary particle";
266  }
267 #endif
268 
269  currentHit = new BscG4Hit;
274 
275  currentHit->setPabs(float(preStepPoint->GetMomentum().mag() * invgev));
277  currentHit->setParticleType(theTrack->GetDefinition()->GetPDGEncoding());
278 
279  float ThetaAtEntry = hitPointLocal.theta() * invdeg;
280  float PhiAtEntry = hitPointLocal.phi() * invdeg;
281 
282  currentHit->setThetaAtEntry(ThetaAtEntry);
283  currentHit->setPhiAtEntry(PhiAtEntry);
284 
288 
289  currentHit->setParentId(theTrack->GetParentID());
290  currentHit->setProcessId(theEnumerator->processId(theTrack->GetCreatorProcess()));
291 
292  currentHit->setVertexPosition(theTrack->GetVertexPosition());
293 
294  updateHit();
296 }
297 
300 
301 #ifdef EDM_ML_DEBUG
302  edm::LogVerbatim("TimingSim") << "updateHit: " << GetName() << " add eloss(GeV) " << edeposit
303  << "CurrentHit=" << currentHit << ", PostStepPoint= " << postStepPoint->GetPosition();
304 #endif
305 
306  // buffer for next steps:
307  tsID = tSliceID;
308  primID = primaryID;
310 }
311 
312 void TimingSD::setToLocal(const G4StepPoint* stepPoint, const G4ThreeVector& globalPoint, G4ThreeVector& localPoint) {
313  const G4VTouchable* touch = stepPoint->GetTouchable();
314  localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(globalPoint);
315 }
316 
317 void TimingSD::EndOfEvent(G4HCofThisEvent*) {
318  int nhits = theHC->entries();
319  if (0 == nhits) {
320  return;
321  }
322  // here we loop over transient hits and make them persistent
323  for (int j = 0; j < nhits; ++j) {
324  BscG4Hit* aHit = (*theHC)[j];
325  Local3DPoint locEntryPoint = ConvertToLocal3DPoint(aHit->getEntryLocalP());
326  Local3DPoint locExitPoint = ConvertToLocal3DPoint(aHit->getExitLocalP());
327 
328 #ifdef EDM_ML_DEBUG
329  edm::LogVerbatim("TimingSim") << "TimingSD: Hit for storage \n"
330  << *aHit << "\n Entry point: " << locEntryPoint << "\n Exit point: " << locExitPoint;
331 #endif
332 
333  slave->processHits(PSimHit(locEntryPoint,
334  locExitPoint,
335  aHit->getPabs(),
336  aHit->getTof(),
337  aHit->getEnergyLoss(),
338  aHit->getParticleType(),
339  aHit->getUnitID(),
340  aHit->getTrackID(),
341  aHit->getThetaAtEntry(),
342  aHit->getPhiAtEntry(),
343  aHit->getProcessId()));
344  }
345 }
346 
348 #ifdef EDM_ML_DEBUG
349  edm::LogVerbatim("TimingSim") << "TimingSD: Collection " << theHC->GetName();
350 #endif
351  theHC->PrintAllHits();
352 }
353 
355  if (slave->name() == hname) {
356  cc = slave->hits();
357  }
358 }
359 
361 #ifdef EDM_ML_DEBUG
362  edm::LogVerbatim("TimingSim") << " Dispatched BeginOfEvent for " << GetName();
363 #endif
364  clearHits();
365 }
366 
uint32_t previousUnitID
Definition: TimingSD.h:83
Log< level::Info, true > LogVerbatim
void setTof(float e)
Definition: BscG4Hit.h:73
void setEntryLocalP(const G4ThreeVector &xyz)
Definition: BscG4Hit.h:33
std::string name() const
static const TGPicture * info(bool iBackgroundIsBlack)
void setEntry(const G4ThreeVector &xyz)
Definition: BSCG4Hit.cc:98
void setTimeFactor(double)
Definition: TimingSD.cc:79
void EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: TimingSD.cc:317
void clearHits() override
Definition: TimingSD.cc:367
double timeFactor
Definition: TimingSD.h:97
void setIncidentEnergy(float e)
Definition: BscG4Hit.h:51
int hcID
Definition: TimingSD.h:86
float edepositHAD
Definition: TimingSD.h:105
uint32_t getUnitID() const
Definition: BscG4Hit.h:56
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
G4ThreeVector hitPointLocal
Definition: TimingSD.h:93
void fillHits(edm::PSimHitContainer &, const std::string &) override
Definition: TimingSD.cc:354
int getParticleType() const
Definition: BscG4Hit.h:70
void setTrackID(int id)
Definition: BscG4Hit.h:54
void setParentId(int p)
Definition: BscG4Hit.h:99
int primID
Definition: TimingSD.h:85
virtual uint32_t setDetUnitId(const G4Step *step)=0
virtual bool checkHit(const G4Step *, BscG4Hit *)
Definition: TimingSD.cc:214
G4ThreeVector hitPointLocalExit
Definition: TimingSD.h:94
void createNewHit(const G4Step *)
Definition: TimingSD.cc:252
int getProcessId() const
Definition: BscG4Hit.h:94
BscG4Hit * currentHit
Definition: TimingSD.h:78
float getPhiAtEntry() const
Definition: BscG4Hit.h:78
float edepositEM
Definition: TimingSD.h:105
void setCuts(double eCut, double historyCut)
Definition: TimingSD.cc:89
void setParticleType(int i)
Definition: BscG4Hit.h:75
float getTof() const
Definition: BscG4Hit.h:68
void setToLocal(const G4StepPoint *stepPoint, const G4ThreeVector &globalPoint, G4ThreeVector &localPoint)
Definition: TimingSD.cc:312
void setThetaAtEntry(float t)
Definition: BscG4Hit.h:80
TimingSD(const std::string &, const SensitiveDetectorCatalog &, const SimTrackManager *)
Definition: TimingSD.cc:37
void storeHit(BscG4Hit *)
Definition: TimingSD.cc:241
float getPabs() const
Definition: BscG4Hit.h:67
void setPabs(float e)
Definition: BscG4Hit.h:72
std::vector< PSimHit > & hits()
void getStepInfo(const G4Step *)
Definition: TimingSD.cc:113
static const float invgev
Definition: TimingSD.cc:32
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
Definition: TimingSD.cc:360
void setUnitID(uint32_t id)
Definition: BscG4Hit.h:57
float tof
Definition: TimingSD.h:103
virtual void Initialize()
bool ProcessHits(G4Step *, G4TouchableHistory *) override
Definition: TimingSD.cc:102
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
void setPhiAtEntry(float f)
Definition: BscG4Hit.h:81
int tsID
Definition: TimingSD.h:87
T sqrt(T t)
Definition: SSEVec.h:19
int getTrackID() const
Definition: BscG4Hit.h:53
float edeposit
Definition: TimingSD.h:104
TrackingSlaveSD * slave
Definition: TimingSD.h:72
G4ProcessTypeEnumerator * theEnumerator
Definition: TimingSD.h:73
static const double invdeg
Definition: TimingSD.cc:34
double f[11][100]
bool hitExists(const G4Step *)
Definition: TimingSD.cc:182
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
const G4Track * theTrack
Definition: TimingSD.h:79
double energyHistoryCut
Definition: TimingSD.h:100
void PrintAll() override
Definition: TimingSD.cc:347
const G4StepPoint * preStepPoint
Definition: TimingSD.h:80
~TimingSD() override
Definition: TimingSD.cc:59
int tSliceID
Definition: TimingSD.h:89
void setTimeSlice(double d)
Definition: BscG4Hit.h:60
void addEnergyDeposit(float em, float hd)
Definition: BSCG4Hit.cc:113
double energyCut
Definition: TimingSD.h:99
float getThetaAtEntry() const
Definition: BscG4Hit.h:77
uint32_t unitID
Definition: TimingSD.h:83
void setExitLocalP(const G4ThreeVector &xyz)
Definition: BscG4Hit.h:36
double tSlice
Definition: TimingSD.h:96
void setVertexPosition(const G4ThreeVector &)
Definition: BSCG4Hit.cc:125
static const double invns
Definition: TimingSD.cc:33
void updateHit()
Definition: TimingSD.cc:298
float getEnergyLoss() const
Definition: BscG4Hit.h:69
int getTimeSliceID() const
Definition: BscG4Hit.h:61
const G4StepPoint * postStepPoint
Definition: TimingSD.h:81
BscG4HitCollection * theHC
Definition: TimingSD.h:76
static bool isGammaElectronPositron(int pdgCode)
const G4ThreeVector & getExitLocalP() const
Definition: BscG4Hit.h:35
virtual bool processHits(const PSimHit &)
std::vector< PSimHit > PSimHitContainer
double incidentEnergy
Definition: TimingSD.h:102
unsigned int processId(const G4VProcess *p) const
G4THitsCollection< BscG4Hit > BscG4HitCollection
Log< level::Warning, false > LogWarning
void setProcessId(int p)
Definition: BscG4Hit.h:100
const G4ThreeVector & getEntryLocalP() const
Definition: BscG4Hit.h:32
void Initialize(G4HCofThisEvent *HCE) override
Definition: TimingSD.cc:64
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
G4ThreeVector hitPointExit
Definition: TimingSD.h:92
G4ThreeVector hitPoint
Definition: TimingSD.h:91
int primaryID
Definition: TimingSD.h:88
#define LogDebug(id)