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