CMS 3D CMS Logo

FP420SD.cc
Go to the documentation of this file.
1 // File: FP420SD.cc
3 // Date: 02.2006
4 // Description: Sensitive Detector class for FP420
5 // Modifications:
7 
12 
16 
21 
25 
28 
29 #include "G4Track.hh"
30 #include "G4SDManager.hh"
31 #include "G4VProcess.hh"
32 #include "G4EventManager.hh"
33 #include "G4Step.hh"
34 #include "G4ParticleTable.hh"
35 
36 #include "G4SystemOfUnits.hh"
37 
38 #include <string>
39 #include <vector>
40 #include <iostream>
41 
43 
44 
45 //#define debug
46 //-------------------------------------------------------------------
48  const SensitiveDetectorCatalog & clg,
49  edm::ParameterSet const & p, const SimTrackManager* manager) :
50  SensitiveTkDetector(name, cpv, clg, p), numberingScheme(nullptr),
51  hcID(-1), theHC(nullptr), theManager(manager), currentHit(nullptr), theTrack(nullptr),
52  currentPV(nullptr), unitID(0), previousUnitID(0), preStepPoint(nullptr),
53  postStepPoint(nullptr), eventno(0){
54 
55  //Parameters
57  int verbn = m_p.getUntrackedParameter<int>("Verbosity");
58  //int verbn = 1;
59 
60  SetVerboseLevel(verbn);
61 
62  slave = new TrackingSlaveSD(name);
63 
64  if (name == "FP420SI") {
65  if (verbn > 0) {
66  edm::LogInfo("FP420Sim") << "name = FP420SI and new FP420NumberingSchem";
67  }
69  } else {
70  edm::LogWarning("FP420Sim") << "FP420SD: ReadoutName not supported\n";
71  }
72 }
73 
75  delete slave;
76  delete numberingScheme;
77 }
78 
79 double FP420SD::getEnergyDeposit(G4Step* aStep) {
80  return aStep->GetTotalEnergyDeposit();
81 }
82 
83 void FP420SD::Initialize(G4HCofThisEvent * HCE) {
84 #ifdef debug
85  LogDebug("FP420Sim") << "FP420SD : Initialize called for " << name << std::endl;
86 #endif
87 
88  theHC = new FP420G4HitCollection(GetName(), collectionName[0]);
89  if (hcID<0)
90  hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
91  HCE->AddHitsCollection(hcID, theHC);
92 
93  tsID = -2;
94  // primID = -2;
95  primID = 0;
96 
98 }
99 
100 
101 bool FP420SD::ProcessHits(G4Step * aStep, G4TouchableHistory * ) {
102 
103  if (aStep == nullptr) {
104  return true;
105  } else {
106  GetStepInfo(aStep);
107  // LogDebug("FP420Sim") << edeposit <<std::endl;
108 
109  //AZ
110 #ifdef debug
111  LogDebug("FP420Sim") << "FP420SD : number of hits = " << theHC->entries() << std::endl;
112 #endif
113 
114  if (HitExists() == false && edeposit>0. && theHC->entries()< 200 ){
115  CreateNewHit();
116  return true;
117  }
118  }
119  return true;
120 }
121 
122 void FP420SD::GetStepInfo(G4Step* aStep) {
123 
124  preStepPoint = aStep->GetPreStepPoint();
125  postStepPoint= aStep->GetPostStepPoint();
126  theTrack = aStep->GetTrack();
127  hitPoint = preStepPoint->GetPosition();
128  currentPV = preStepPoint->GetPhysicalVolume();
129  hitPointExit = postStepPoint->GetPosition();
130 
131  hitPointLocal = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
132  hitPointLocalExit = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPointExit);
133 
134 
135  G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
136  if (particleCode == emPDG ||
137  particleCode == epPDG ||
138  particleCode == gammaPDG ) {
139  edepositEM = getEnergyDeposit(aStep); edepositHAD = 0.;
140  } else {
141  edepositEM = 0.; edepositHAD = getEnergyDeposit(aStep);
142  }
143  edeposit = aStep->GetTotalEnergyDeposit();
144  tSlice = (postStepPoint->GetGlobalTime() )/nanosecond;
145  tSliceID = (int) tSlice;
146  unitID = setDetUnitId(aStep);
147 #ifdef debug
148  LogDebug("FP420Sim") << "unitID=" << unitID <<std::endl;
149 #endif
150  primaryID = theTrack->GetTrackID();
151  // Position = hitPoint;
152  Pabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
153  //Tof = 1400. + aStep->GetPostStepPoint()->GetGlobalTime()/nanosecond;
154  Tof = aStep->GetPostStepPoint()->GetGlobalTime()/nanosecond;
155  Eloss = aStep->GetTotalEnergyDeposit()/GeV;
156  ParticleType = theTrack->GetDefinition()->GetPDGEncoding();
157  ThetaAtEntry = aStep->GetPreStepPoint()->GetPosition().theta()/deg;
158  PhiAtEntry = aStep->GetPreStepPoint()->GetPosition().phi()/deg;
159 
160  ParentId = theTrack->GetParentID();
161  Vx = theTrack->GetVertexPosition().x();
162  Vy = theTrack->GetVertexPosition().y();
163  Vz = theTrack->GetVertexPosition().z();
164  X = hitPoint.x();
165  Y = hitPoint.y();
166  Z = hitPoint.z();
167 }
168 
169 uint32_t FP420SD::setDetUnitId(const G4Step * aStep) {
170 
171  return (numberingScheme == nullptr ? 0 : numberingScheme->getUnitID(aStep));
172 }
173 
175  if (primaryID<1) {
176  edm::LogWarning("FP420Sim") << "***** FP420SD error: primaryID = "
177  << primaryID
178  << " maybe detector name changed";
179  }
180 
181  // Update if in the same detector, time-slice and for same track
182  // if (primaryID == primID && tSliceID == tsID && unitID==previousUnitID) {
183  if (tSliceID == tsID && unitID==previousUnitID) {
184  //AZ:
185  UpdateHit();
186  return true;
187  }
188  // Reset entry point for new primary
189  if (primaryID != primID)
191 
192  //look in the HitContainer whether a hit with the same primID, unitID,
193  //tSliceID already exists:
194 
195  G4bool found = false;
196 
197  // LogDebug("FP420Sim") << "FP420SD: HCollection= " << theHC->entries() <<std::endl;
198 
199  for (int j=0; j<theHC->entries()&&!found; j++) {
200  FP420G4Hit* aPreviousHit = (*theHC)[j];
201  if (aPreviousHit->getTrackID() == primaryID &&
202  aPreviousHit->getTimeSliceID() == tSliceID &&
203  aPreviousHit->getUnitID() == unitID ) {
204  //AZ:
205  currentHit = aPreviousHit;
206  found = true;
207  }
208  }
209 
210  if (found) {
211  //AZ:
212  UpdateHit();
213  return true;
214  } else {
215  return false;
216  }
217 }
218 
219 
221 
224  incidentEnergy = preStepPoint->GetKineticEnergy();
225 
226 }
227 
228 
230 
231  // if (primID<0) return;
232  if (hit == nullptr ) {
233  edm::LogWarning("FP420Sim") << "FP420SD: hit to be stored is NULL !!";
234  return;
235  }
236 
237  theHC->insert( hit );
238 }
239 
240 
242 
243 #ifdef debug
244  // << " MVid = " << currentPV->GetMother()->GetCopyNo()
245  LogDebug("FP420Sim") << "FP420SD CreateNewHit for"
246  << " PV " << currentPV->GetName()
247  << " PVid = " << currentPV->GetCopyNo()
248  << " Unit " << unitID <<std::endl;
249  LogDebug("FP420Sim") << " primary " << primaryID
250  << " time slice " << tSliceID
251  << " For Track " << theTrack->GetTrackID()
252  << " which is a " << theTrack->GetDefinition()->GetParticleName();
253 
254  if (theTrack->GetTrackID()==1) {
255  LogDebug("FP420Sim") << " of energy " << theTrack->GetTotalEnergy();
256  } else {
257  LogDebug("FP420Sim") << " daughter of part. " << theTrack->GetParentID();
258  }
259 
260  LogDebug("FP420Sim") << " and created by " ;
261  if (theTrack->GetCreatorProcess()!=NULL)
262  LogDebug("FP420Sim") << theTrack->GetCreatorProcess()->GetProcessName() ;
263  else
264  LogDebug("FP420Sim") << "NO process";
265  LogDebug("FP420Sim") << std::endl;
266 #endif
267 
268 
269  currentHit = new FP420G4Hit;
274 
281 
282 // currentHit->setEntry(entrancePoint);
284 
287 
289  currentHit->setVx(Vx);
290  currentHit->setVy(Vy);
291  currentHit->setVz(Vz);
292 
293  currentHit->setX(X);
294  currentHit->setY(Y);
295  currentHit->setZ(Z);
296  //AZ:12.10.2007
297  // UpdateHit();
298  // buffer for next steps:
299  tsID = tSliceID;
300  primID = primaryID;
302 
304 }
305 
306 
308  //
309  if (Eloss > 0.) {
311 
312 #ifdef debug
313  LogDebug("FP420Sim") << "updateHit: add eloss " << Eloss <<std::endl;
314  LogDebug("FP420Sim") << "CurrentHit=" << currentHit
315  << ", PostStepPoint=" << postStepPoint->GetPosition() << std::endl;
316 #endif
317  //AZ
318  // currentHit->setEnergyLoss(Eloss);
320  }
321 
322 
323 
324  // buffer for next steps:
325  tsID = tSliceID;
326  primID = primaryID;
328 }
329 
330 
331 G4ThreeVector FP420SD::SetToLocal(const G4ThreeVector& global){
332 
333  const G4VTouchable* touch= preStepPoint->GetTouchable();
334  theEntryPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
335  return theEntryPoint;
336 }
337 
338 
339 G4ThreeVector FP420SD::SetToLocalExit(const G4ThreeVector& globalPoint){
340 
341  const G4VTouchable* touch= postStepPoint->GetTouchable();
342  theExitPoint = touch->GetHistory()->GetTopTransform().TransformPoint(globalPoint);
343  return theExitPoint;
344 }
345 
346 
347 void FP420SD::EndOfEvent(G4HCofThisEvent* ) {
348 
349  // here we loop over transient hits and make them persistent
350 
351  // if(theHC->entries() > 100){
352  // LogDebug("FP420Sim") << "FP420SD: warning!!! Number of hits exceed 100 and =" << theHC->entries() << "\n";
353  // }
354  // for (int j=0; j<theHC->entries() && j<100; j++) {
355  int nhitsHPS240=0;
356  int nhitsFP420=0;
357  for (int j=0; j<theHC->entries(); j++) {
358  FP420G4Hit* aHit = (*theHC)[j];
359  if((fabs(aHit->getTof())> 780. && fabs(aHit->getTof())< 840.)) ++nhitsHPS240;
360  if((fabs(aHit->getTof())>1380. && fabs(aHit->getTof())<1450.)) ++nhitsFP420;
361  // if(fabs(aHit->getTof()) < 1700.) {
362  if((fabs(aHit->getTof())>780. && fabs(aHit->getTof())<840. && nhitsHPS240<200.) || ( fabs(aHit->getTof())>1380. && fabs(aHit->getTof())<1450. && nhitsFP420<200.)) {
363 
364 
365 #ifdef ddebug
366  // LogDebug("FP420SD") << " FP420Hit " << j << " " << *aHit << std::endl;
367  LogDebug("FP420Sim") << "hit number" << j << "unit ID = "<<aHit->getUnitID()<< "\n";
368  LogDebug("FP420Sim") << "entry z " << aHit->getEntry().z()<< "\n";
369  LogDebug("FP420Sim") << "entr theta " << aHit->getThetaAtEntry()<< "\n";
370 #endif
371 
372  // Local3DPoint locExitPoint(0,0,0);
373  // Local3DPoint locEntryPoint(aHit->getEntry().x(),
374  // aHit->getEntry().y(),
375  // aHit->getEntry().z());
376  Local3DPoint locExitPoint(aHit->getExitLocalP().x(),
377  aHit->getExitLocalP().y(),
378  aHit->getExitLocalP().z());
379  Local3DPoint locEntryPoint(aHit->getEntryLocalP().x(),
380  aHit->getEntryLocalP().y(),
381  aHit->getEntryLocalP().z());
382  // implicit conversion (slicing) to PSimHit!!!
383  slave->processHits(PSimHit(locEntryPoint,locExitPoint,//entryPoint(), exitPoint() Local3DPoint
384  aHit->getPabs(),// pabs() float
385  aHit->getTof(), // tof() float
386  aHit->getEnergyLoss(), // energyLoss() float
387  aHit->getParticleType(),// particleType() int
388  aHit->getUnitID(), // detUnitId() unsigned int
389  aHit->getTrackID(),// trackId() unsigned int
390  aHit->getThetaAtEntry(),// thetaAtEntry() float
391  aHit->getPhiAtEntry())); // phiAtEntry() float
392 
393  //PSimHit( const Local3DPoint& entry, const Local3DPoint& exit,
394  // float pabs, float tof, float eloss, int particleType,
395  // unsigned int detId, unsigned int trackId,
396  // float theta, float phi, unsigned short processType=0) :
397 
398  // LocalVector direction = hit.exitPoint() - hit.entryPoint();
399  //hit.energyLoss
400 
401 
402  /*
403  aHit->getEM(), -
404  aHit->getHadr(), -
405  aHit->getIncidentEnergy(), -
406  aHit->getTimeSlice(), -
407  aHit->getEntry(), -
408  aHit->getParentId(),
409  aHit->getEntryLocalP(), -
410  aHit->getExitLocalP(), -
411  aHit->getX(), -
412  aHit->getY(), -
413  aHit->getZ(), -
414  aHit->getVx(), -
415  aHit->getVy(), -
416  aHit->getVz())); -
417  */
418  }//if Tof<1600. if nhits<100
419  } // for loop on hits
420 
421  Summarize();
422 }
423 
424 
426 }
427 
428 
430 }
431 
432 
434 }
435 
436 
438  LogDebug("FP420Sim") << "FP420SD: Collection " << theHC->GetName() << "\n";
439  theHC->PrintAllHits();
440 }
441 
442 
443 //void FP420SD::SetNumberingScheme(FP420NumberingScheme* scheme){
444 //
445 // if (numberingScheme)
446 // delete numberingScheme;
447 // numberingScheme = scheme;
448 //
449 //}
450 
452  if (slave->name() == hname) { cc=slave->hits(); }
453 }
454 
456  LogDebug("ForwardSim") << " Dispatched BeginOfEvent for " << GetName()
457  << " !" ;
458  clearHits();
459  eventno = (*i)()->GetEventID();
460 }
461 
463 
464  G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
465  G4String particleName;
466  emPDG = theParticleTable->FindParticle(particleName="e-")->GetPDGEncoding();
467  epPDG = theParticleTable->FindParticle(particleName="e+")->GetPDGEncoding();
468  gammaPDG = theParticleTable->FindParticle(particleName="gamma")->GetPDGEncoding();
469 
470 }
471 
472 void FP420SD::update (const ::EndOfEvent*) {
473 }
474 
476  //AZ:
477  slave->Initialize();
478 }
479 
#define LogDebug(id)
int eventno
Definition: FP420SD.h:133
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
G4ThreeVector theExitPoint
Definition: FP420SD.h:90
G4int emPDG
Definition: FP420SD.h:138
void setY(float t)
Definition: FP420G4Hit.cc:169
const double GeV
Definition: MathUtil.h:16
G4ThreeVector SetToLocalExit(const G4ThreeVector &globalPoint)
Definition: FP420SD.cc:339
void setPhiAtEntry(float f)
Definition: FP420G4Hit.cc:163
float getTof() const
Definition: FP420G4Hit.cc:149
void DrawAll() override
Definition: FP420SD.cc:433
float ThetaAtEntry
Definition: FP420SD.h:122
void setVz(float p)
Definition: FP420G4Hit.cc:184
void setEntry(const G4ThreeVector &xyz)
Definition: FP420G4Hit.cc:116
virtual double getEnergyDeposit(G4Step *step)
Definition: FP420SD.cc:79
G4ThreeVector hitPointExit
Definition: FP420SD.h:114
std::string name() const
G4ThreeVector getExitLocalP() const
Definition: FP420G4Hit.cc:121
void setZ(float t)
Definition: FP420G4Hit.cc:172
float Pabs
Definition: FP420SD.h:117
void setPabs(float e)
Definition: FP420G4Hit.cc:153
TrackingSlaveSD * slave
Definition: FP420SD.h:85
#define NULL
Definition: scimark2.h:8
G4StepPoint * preStepPoint
Definition: FP420SD.h:109
void ResetForNewPrimary()
Definition: FP420SD.cc:220
float getPabs() const
Definition: FP420G4Hit.cc:148
#define nullptr
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
void setVx(float p)
Definition: FP420G4Hit.cc:178
G4int hcID
Definition: FP420SD.h:94
void setTrackID(int i)
Definition: FP420G4Hit.cc:134
unsigned int getTrackID() const
Definition: FP420G4Hit.cc:133
void addEnergyDeposit(double em, double hd)
Definition: FP420G4Hit.cc:143
G4int tsID
Definition: FP420SD.h:98
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
unsigned int primaryID
Definition: FP420SD.h:105
float Eloss
Definition: FP420SD.h:119
unsigned int primID
Definition: FP420SD.h:105
float getPhiAtEntry() const
Definition: FP420G4Hit.cc:160
G4ThreeVector exitPoint
Definition: FP420SD.h:88
float edeposit
Definition: FP420SD.h:111
unsigned int getUnitID() const
Definition: FP420G4Hit.cc:136
void setEnergyLoss(float e)
Definition: FP420G4Hit.cc:156
std::vector< PSimHit > & hits()
void clearHits() override
Definition: FP420SD.cc:475
std::string const collectionName[nCollections]
Definition: Collections.h:47
G4ThreeVector theEntryPoint
Definition: FP420SD.h:89
void setTof(float e)
Definition: FP420G4Hit.cc:154
virtual void Initialize()
float edepositHAD
Definition: FP420SD.h:137
float Vz
Definition: FP420SD.h:126
void GetStepInfo(G4Step *aStep)
Definition: FP420SD.cc:122
void fillHits(edm::PSimHitContainer &, const std::string &) override
Definition: FP420SD.cc:451
G4int gammaPDG
Definition: FP420SD.h:140
FP420G4HitCollection * theHC
Definition: FP420SD.h:95
void CreateNewHit()
Definition: FP420SD.cc:241
G4Track * theTrack
Definition: FP420SD.h:100
void Summarize()
Definition: FP420SD.cc:425
G4int tSliceID
Definition: FP420SD.h:104
G4ThreeVector getEntryLocalP() const
Definition: FP420G4Hit.cc:118
G4bool HitExists()
Definition: FP420SD.cc:174
void update(const BeginOfRun *) override
This routine will be called when the appropriate signal arrives.
Definition: FP420SD.cc:462
void setTimeSlice(double d)
Definition: FP420G4Hit.cc:140
virtual unsigned int getUnitID(const G4Step *aStep) const
float getEnergyLoss() const
Definition: FP420G4Hit.cc:150
void setEntryLocalP(const G4ThreeVector &xyz1)
Definition: FP420G4Hit.cc:119
uint32_t setDetUnitId(const G4Step *) override
Definition: FP420SD.cc:169
G4THitsCollection< FP420G4Hit > FP420G4HitCollection
void setIncidentEnergy(double e)
Definition: FP420G4Hit.cc:131
uint32_t previousUnitID
Definition: FP420SD.h:103
int getTimeSliceID() const
Definition: FP420G4Hit.cc:141
G4VPhysicalVolume * currentPV
Definition: FP420SD.h:101
void Initialize(G4HCofThisEvent *HCE) override
Definition: FP420SD.cc:83
FP420SD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: FP420SD.cc:47
void setX(float t)
Definition: FP420G4Hit.cc:166
float Tof
Definition: FP420SD.h:118
void setExitLocalP(const G4ThreeVector &xyz1)
Definition: FP420G4Hit.cc:122
void clear() override
Definition: FP420SD.cc:429
int getParticleType() const
Definition: FP420G4Hit.cc:151
uint32_t unitID
Definition: FP420SD.h:103
float PhiAtEntry
Definition: FP420SD.h:123
float edepositEM
Definition: FP420SD.h:137
short ParticleType
Definition: FP420SD.h:120
void setParticleType(short i)
Definition: FP420G4Hit.cc:157
void setThetaAtEntry(float t)
Definition: FP420G4Hit.cc:162
G4ThreeVector hitPointLocalExit
Definition: FP420SD.h:116
FP420NumberingScheme * numberingScheme
Definition: FP420SD.h:86
G4double tSlice
Definition: FP420SD.h:107
bool ProcessHits(G4Step *, G4TouchableHistory *) override
Definition: FP420SD.cc:101
void PrintAll() override
Definition: FP420SD.cc:437
G4ThreeVector SetToLocal(const G4ThreeVector &global)
Definition: FP420SD.cc:331
float incidentEnergy
Definition: FP420SD.h:92
void setUnitID(unsigned int i)
Definition: FP420G4Hit.cc:137
int ParentId
Definition: FP420SD.h:125
G4ThreeVector entrancePoint
Definition: FP420SD.h:88
float Z
Definition: FP420SD.h:127
float Vy
Definition: FP420SD.h:126
virtual bool processHits(const PSimHit &)
void UpdateHit()
Definition: FP420SD.cc:307
std::vector< PSimHit > PSimHitContainer
G4ThreeVector getEntry() const
Definition: FP420G4Hit.cc:115
~FP420SD() override
Definition: FP420SD.cc:74
G4int epPDG
Definition: FP420SD.h:139
G4StepPoint * postStepPoint
Definition: FP420SD.h:110
void setParentId(int p)
Definition: FP420G4Hit.cc:175
float X
Definition: FP420SD.h:127
FP420G4Hit * currentHit
Definition: FP420SD.h:99
G4ThreeVector hitPointLocal
Definition: FP420SD.h:115
G4ThreeVector hitPoint
Definition: FP420SD.h:113
void addEnergyLoss(float e)
Definition: FP420G4Hit.cc:155
float getThetaAtEntry() const
Definition: FP420G4Hit.cc:159
float Y
Definition: FP420SD.h:127
void StoreHit(FP420G4Hit *)
Definition: FP420SD.cc:229
void setVy(float p)
Definition: FP420G4Hit.cc:181
float Vx
Definition: FP420SD.h:126
void EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: FP420SD.cc:347