forked from cms-sw/cmssw
-
Notifications
You must be signed in to change notification settings - Fork 0
/
FP420SD.cc
445 lines (377 loc) · 14.5 KB
/
FP420SD.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
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
///////////////////////////////////////////////////////////////////////////////
// File: FP420SD.cc
// Date: 02.2006
// Description: Sensitive Detector class for FP420
// Modifications:
///////////////////////////////////////////////////////////////////////////////
#include "SimG4Core/SensitiveDetector/interface/FrameRotation.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 "SimDataFormats/TrackingHit/interface/UpdatablePSimHit.h"
#include "DataFormats/GeometryVector/interface/LocalPoint.h"
#include "SimG4CMS/FP420/interface/FP420SD.h"
#include "SimG4CMS/FP420/interface/FP420G4Hit.h"
#include "SimG4CMS/FP420/interface/FP420G4HitCollection.h"
#include "SimG4CMS/FP420/interface/FP420NumberingScheme.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "SimG4Core/Notification/interface/TrackInformation.h"
#include "SimG4Core/Notification/interface/SimTrackManager.h"
#include "G4Track.hh"
#include "G4SDManager.hh"
#include "G4VProcess.hh"
#include "G4EventManager.hh"
#include "G4Step.hh"
#include "G4ParticleTable.hh"
#include "G4SystemOfUnits.hh"
#include <string>
#include <vector>
#include <iostream>
#include "FWCore/MessageLogger/interface/MessageLogger.h"
//#define debug
//-------------------------------------------------------------------
FP420SD::FP420SD(const std::string& name,
const edm::EventSetup& es,
const SensitiveDetectorCatalog& clg,
edm::ParameterSet const& p,
const SimTrackManager* manager)
: SensitiveTkDetector(name, es, clg, p),
numberingScheme(nullptr),
hcID(-1),
theHC(nullptr),
theManager(manager),
currentHit(nullptr),
theTrack(nullptr),
currentPV(nullptr),
unitID(0),
previousUnitID(0),
preStepPoint(nullptr),
postStepPoint(nullptr),
eventno(0) {
//Parameters
edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("FP420SD");
int verbn = m_p.getUntrackedParameter<int>("Verbosity");
//int verbn = 1;
SetVerboseLevel(verbn);
slave = new TrackingSlaveSD(name);
if (name == "FP420SI") {
if (verbn > 0) {
edm::LogInfo("FP420Sim") << "name = FP420SI and new FP420NumberingSchem";
}
numberingScheme = new FP420NumberingScheme();
} else {
edm::LogWarning("FP420Sim") << "FP420SD: ReadoutName not supported\n";
}
}
FP420SD::~FP420SD() {
delete slave;
delete numberingScheme;
}
double FP420SD::getEnergyDeposit(G4Step* aStep) { return aStep->GetTotalEnergyDeposit(); }
void FP420SD::Initialize(G4HCofThisEvent* HCE) {
#ifdef debug
LogDebug("FP420Sim") << "FP420SD : Initialize called for " << name << std::endl;
#endif
theHC = new FP420G4HitCollection(GetName(), collectionName[0]);
if (hcID < 0)
hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
HCE->AddHitsCollection(hcID, theHC);
tsID = -2;
// primID = -2;
primID = 0;
//// slave->Initialize();
}
bool FP420SD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
if (aStep == nullptr) {
return true;
} else {
GetStepInfo(aStep);
// LogDebug("FP420Sim") << edeposit <<std::endl;
//AZ
#ifdef debug
LogDebug("FP420Sim") << "FP420SD : number of hits = " << theHC->entries() << std::endl;
#endif
if (HitExists() == false && edeposit > 0. && theHC->entries() < 200) {
CreateNewHit();
return true;
}
}
return true;
}
void FP420SD::GetStepInfo(G4Step* aStep) {
preStepPoint = aStep->GetPreStepPoint();
postStepPoint = aStep->GetPostStepPoint();
theTrack = aStep->GetTrack();
hitPoint = preStepPoint->GetPosition();
currentPV = preStepPoint->GetPhysicalVolume();
hitPointExit = postStepPoint->GetPosition();
hitPointLocal = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
hitPointLocalExit = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPointExit);
G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
if (particleCode == emPDG || particleCode == epPDG || particleCode == gammaPDG) {
edepositEM = getEnergyDeposit(aStep);
edepositHAD = 0.;
} else {
edepositEM = 0.;
edepositHAD = getEnergyDeposit(aStep);
}
edeposit = aStep->GetTotalEnergyDeposit();
tSlice = (postStepPoint->GetGlobalTime()) / nanosecond;
tSliceID = (int)tSlice;
unitID = setDetUnitId(aStep);
#ifdef debug
LogDebug("FP420Sim") << "unitID=" << unitID << std::endl;
#endif
primaryID = theTrack->GetTrackID();
// Position = hitPoint;
Pabs = aStep->GetPreStepPoint()->GetMomentum().mag() / GeV;
//Tof = 1400. + aStep->GetPostStepPoint()->GetGlobalTime()/nanosecond;
Tof = aStep->GetPostStepPoint()->GetGlobalTime() / nanosecond;
Eloss = aStep->GetTotalEnergyDeposit() / GeV;
ParticleType = theTrack->GetDefinition()->GetPDGEncoding();
ThetaAtEntry = aStep->GetPreStepPoint()->GetPosition().theta() / deg;
PhiAtEntry = aStep->GetPreStepPoint()->GetPosition().phi() / deg;
ParentId = theTrack->GetParentID();
Vx = theTrack->GetVertexPosition().x();
Vy = theTrack->GetVertexPosition().y();
Vz = theTrack->GetVertexPosition().z();
X = hitPoint.x();
Y = hitPoint.y();
Z = hitPoint.z();
}
uint32_t FP420SD::setDetUnitId(const G4Step* aStep) {
return (numberingScheme == nullptr ? 0 : numberingScheme->getUnitID(aStep));
}
G4bool FP420SD::HitExists() {
if (primaryID < 1) {
edm::LogWarning("FP420Sim") << "***** FP420SD error: primaryID = " << primaryID << " maybe detector name changed";
}
// Update if in the same detector, time-slice and for same track
// if (primaryID == primID && tSliceID == tsID && unitID==previousUnitID) {
if (tSliceID == tsID && unitID == previousUnitID) {
//AZ:
UpdateHit();
return true;
}
// Reset entry point for new primary
if (primaryID != primID)
ResetForNewPrimary();
//look in the HitContainer whether a hit with the same primID, unitID,
//tSliceID already exists:
G4bool found = false;
// LogDebug("FP420Sim") << "FP420SD: HCollection= " << theHC->entries() <<std::endl;
int nhits = theHC->entries();
for (int j = 0; j < nhits && !found; j++) {
FP420G4Hit* aPreviousHit = (*theHC)[j];
if (aPreviousHit->getTrackID() == primaryID && aPreviousHit->getTimeSliceID() == tSliceID &&
aPreviousHit->getUnitID() == unitID) {
//AZ:
currentHit = aPreviousHit;
found = true;
}
}
if (found) {
//AZ:
UpdateHit();
return true;
} else {
return false;
}
}
void FP420SD::ResetForNewPrimary() {
entrancePoint = SetToLocal(hitPoint);
exitPoint = SetToLocalExit(hitPointExit);
incidentEnergy = preStepPoint->GetKineticEnergy();
}
void FP420SD::StoreHit(FP420G4Hit* hit) {
// if (primID<0) return;
if (hit == nullptr) {
edm::LogWarning("FP420Sim") << "FP420SD: hit to be stored is NULL !!";
return;
}
theHC->insert(hit);
}
void FP420SD::CreateNewHit() {
#ifdef debug
// << " MVid = " << currentPV->GetMother()->GetCopyNo()
LogDebug("FP420Sim") << "FP420SD CreateNewHit for"
<< " PV " << currentPV->GetName() << " PVid = " << currentPV->GetCopyNo() << " Unit " << unitID
<< std::endl;
LogDebug("FP420Sim") << " primary " << primaryID << " time slice " << tSliceID << " For Track "
<< theTrack->GetTrackID() << " which is a " << theTrack->GetDefinition()->GetParticleName();
if (theTrack->GetTrackID() == 1) {
LogDebug("FP420Sim") << " of energy " << theTrack->GetTotalEnergy();
} else {
LogDebug("FP420Sim") << " daughter of part. " << theTrack->GetParentID();
}
LogDebug("FP420Sim") << " and created by ";
if (theTrack->GetCreatorProcess() != NULL)
LogDebug("FP420Sim") << theTrack->GetCreatorProcess()->GetProcessName();
else
LogDebug("FP420Sim") << "NO process";
LogDebug("FP420Sim") << std::endl;
#endif
currentHit = new FP420G4Hit;
currentHit->setTrackID(primaryID);
currentHit->setTimeSlice(tSlice);
currentHit->setUnitID(unitID);
currentHit->setIncidentEnergy(incidentEnergy);
currentHit->setPabs(Pabs);
currentHit->setTof(Tof);
currentHit->setEnergyLoss(Eloss);
currentHit->setParticleType(ParticleType);
currentHit->setThetaAtEntry(ThetaAtEntry);
currentHit->setPhiAtEntry(PhiAtEntry);
// currentHit->setEntry(entrancePoint);
currentHit->setEntry(hitPoint);
currentHit->setEntryLocalP(hitPointLocal);
currentHit->setExitLocalP(hitPointLocalExit);
currentHit->setParentId(ParentId);
currentHit->setVx(Vx);
currentHit->setVy(Vy);
currentHit->setVz(Vz);
currentHit->setX(X);
currentHit->setY(Y);
currentHit->setZ(Z);
//AZ:12.10.2007
// UpdateHit();
// buffer for next steps:
tsID = tSliceID;
primID = primaryID;
previousUnitID = unitID;
StoreHit(currentHit);
}
void FP420SD::UpdateHit() {
//
if (Eloss > 0.) {
currentHit->addEnergyDeposit(edepositEM, edepositHAD);
#ifdef debug
LogDebug("FP420Sim") << "updateHit: add eloss " << Eloss << std::endl;
LogDebug("FP420Sim") << "CurrentHit=" << currentHit << ", PostStepPoint=" << postStepPoint->GetPosition()
<< std::endl;
#endif
//AZ
// currentHit->setEnergyLoss(Eloss);
currentHit->addEnergyLoss(Eloss);
}
// buffer for next steps:
tsID = tSliceID;
primID = primaryID;
previousUnitID = unitID;
}
G4ThreeVector FP420SD::SetToLocal(const G4ThreeVector& global) {
const G4VTouchable* touch = preStepPoint->GetTouchable();
theEntryPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
return theEntryPoint;
}
G4ThreeVector FP420SD::SetToLocalExit(const G4ThreeVector& globalPoint) {
const G4VTouchable* touch = postStepPoint->GetTouchable();
theExitPoint = touch->GetHistory()->GetTopTransform().TransformPoint(globalPoint);
return theExitPoint;
}
void FP420SD::EndOfEvent(G4HCofThisEvent*) {
// here we loop over transient hits and make them persistent
// if(theHC->entries() > 100){
// LogDebug("FP420Sim") << "FP420SD: warning!!! Number of hits exceed 100 and =" << theHC->entries() << "\n";
// }
// for (int j=0; j<theHC->entries() && j<100; j++) {
int nhitsHPS240 = 0;
int nhitsFP420 = 0;
int nhits = theHC->entries();
for (int j = 0; j < nhits; j++) {
FP420G4Hit* aHit = (*theHC)[j];
if ((fabs(aHit->getTof()) > 780. && fabs(aHit->getTof()) < 840.))
++nhitsHPS240;
if ((fabs(aHit->getTof()) > 1380. && fabs(aHit->getTof()) < 1450.))
++nhitsFP420;
// if(fabs(aHit->getTof()) < 1700.) {
if ((fabs(aHit->getTof()) > 780. && fabs(aHit->getTof()) < 840. && nhitsHPS240 < 200.) ||
(fabs(aHit->getTof()) > 1380. && fabs(aHit->getTof()) < 1450. && nhitsFP420 < 200.)) {
#ifdef ddebug
// LogDebug("FP420SD") << " FP420Hit " << j << " " << *aHit << std::endl;
LogDebug("FP420Sim") << "hit number" << j << "unit ID = " << aHit->getUnitID() << "\n";
LogDebug("FP420Sim") << "entry z " << aHit->getEntry().z() << "\n";
LogDebug("FP420Sim") << "entr theta " << aHit->getThetaAtEntry() << "\n";
#endif
// Local3DPoint locExitPoint(0,0,0);
// Local3DPoint locEntryPoint(aHit->getEntry().x(),
// aHit->getEntry().y(),
// aHit->getEntry().z());
Local3DPoint locExitPoint(aHit->getExitLocalP().x(), aHit->getExitLocalP().y(), aHit->getExitLocalP().z());
Local3DPoint locEntryPoint(aHit->getEntryLocalP().x(), aHit->getEntryLocalP().y(), aHit->getEntryLocalP().z());
// implicit conversion (slicing) to PSimHit!!!
slave->processHits(PSimHit(locEntryPoint,
locExitPoint, //entryPoint(), exitPoint() Local3DPoint
aHit->getPabs(), // pabs() float
aHit->getTof(), // tof() float
aHit->getEnergyLoss(), // energyLoss() float
aHit->getParticleType(), // particleType() int
aHit->getUnitID(), // detUnitId() unsigned int
aHit->getTrackID(), // trackId() unsigned int
aHit->getThetaAtEntry(), // thetaAtEntry() float
aHit->getPhiAtEntry())); // phiAtEntry() float
//PSimHit( const Local3DPoint& entry, const Local3DPoint& exit,
// float pabs, float tof, float eloss, int particleType,
// unsigned int detId, unsigned int trackId,
// float theta, float phi, unsigned short processType=0) :
// LocalVector direction = hit.exitPoint() - hit.entryPoint();
//hit.energyLoss
/*
aHit->getEM(), -
aHit->getHadr(), -
aHit->getIncidentEnergy(), -
aHit->getTimeSlice(), -
aHit->getEntry(), -
aHit->getParentId(),
aHit->getEntryLocalP(), -
aHit->getExitLocalP(), -
aHit->getX(), -
aHit->getY(), -
aHit->getZ(), -
aHit->getVx(), -
aHit->getVy(), -
aHit->getVz())); -
*/
} //if Tof<1600. if nhits<100
} // for loop on hits
Summarize();
}
void FP420SD::Summarize() {}
void FP420SD::clear() {}
void FP420SD::DrawAll() {}
void FP420SD::PrintAll() {
LogDebug("FP420Sim") << "FP420SD: Collection " << theHC->GetName() << "\n";
theHC->PrintAllHits();
}
//void FP420SD::SetNumberingScheme(FP420NumberingScheme* scheme){
//
// if (numberingScheme)
// delete numberingScheme;
// numberingScheme = scheme;
//
//}
void FP420SD::fillHits(edm::PSimHitContainer& cc, const std::string& hname) {
if (slave->name() == hname) {
cc = slave->hits();
}
}
void FP420SD::update(const BeginOfEvent* i) {
LogDebug("ForwardSim") << " Dispatched BeginOfEvent for " << GetName() << " !";
clearHits();
eventno = (*i)()->GetEventID();
}
void FP420SD::update(const BeginOfRun*) {
G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
G4String particleName;
emPDG = theParticleTable->FindParticle(particleName = "e-")->GetPDGEncoding();
epPDG = theParticleTable->FindParticle(particleName = "e+")->GetPDGEncoding();
gammaPDG = theParticleTable->FindParticle(particleName = "gamma")->GetPDGEncoding();
}
void FP420SD::update(const ::EndOfEvent*) {}
void FP420SD::clearHits() {
//AZ:
slave->Initialize();
}