CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BscSD.cc
Go to the documentation of this file.
1 // File: BscSD.cc
3 // Date: 02.2006
4 // Description: Sensitive Detector class for Bsc
5 // Modifications:
7 
8 //#include "Geometry/Vector/interface/LocalPoint.h"
9 //#include "Geometry/Vector/interface/LocalVector.h"
10 
14 
17 
18 
23 
27 
29 
30 #include "G4Track.hh"
31 #include "G4SDManager.hh"
32 #include "G4VProcess.hh"
33 #include "G4EventManager.hh"
34 #include "G4Step.hh"
35 #include "G4ParticleTable.hh"
36 
37 #include "G4PhysicalConstants.hh"
38 #include "G4SystemOfUnits.hh"
39 
40 #include <string>
41 #include <vector>
42 #include <iostream>
43 
45 
46 #define debug
47 //-------------------------------------------------------------------
50  edm::ParameterSet const & p, const SimTrackManager* manager) :
51  SensitiveTkDetector(name, cpv, clg, p), numberingScheme(0), name(name),
52  hcID(-1), theHC(0), theManager(manager), currentHit(0), theTrack(0),
53  currentPV(0), unitID(0), previousUnitID(0), preStepPoint(0),
54  postStepPoint(0), eventno(0){
55 
56  //Add Bsc Sentitive Detector Name
57  collectionName.insert(name);
58 
59 
60  //Parameters
62  int verbn = m_p.getUntrackedParameter<int>("Verbosity");
63  //int verbn = 1;
64 
65  SetVerboseLevel(verbn);
66  LogDebug("BscSim")
67  << "*******************************************************\n"
68  << "* *\n"
69  << "* Constructing a BscSD with name " << name << "\n"
70  << "* *\n"
71  << "*******************************************************";
72 
73 
74  slave = new TrackingSlaveSD(name);
75 
76  //
77  // attach detectors (LogicalVolumes)
78  //
79  std::vector<std::string> lvNames = clg.logicalNames(name);
80 
81  this->Register();
82 
83  for (std::vector<std::string>::iterator it=lvNames.begin();
84  it !=lvNames.end(); it++) {
85  this->AssignSD(*it);
86  edm::LogInfo("BscSim") << "BscSD : Assigns SD to LV " << (*it);
87  }
88 
89  if (name == "BSCHits") {
90  if (verbn > 0) {
91  edm::LogInfo("BscSim") << "name = BSCHits and new BscNumberingSchem";
92  }
94  } else {
95  edm::LogWarning("BscSim") << "BscSD: ReadoutName "<<name<<" not supported";
96  }
97 
98  edm::LogInfo("BscSim") << "BscSD: Instantiation completed";
99 }
100 
101 
102 
103 
105  //AZ:
106  if (slave) delete slave;
107 
108  if (numberingScheme)
109  delete numberingScheme;
110 
111 }
112 
113 double BscSD::getEnergyDeposit(G4Step* aStep) {
114  return aStep->GetTotalEnergyDeposit();
115 }
116 
117 void BscSD::Initialize(G4HCofThisEvent * HCE) {
118 #ifdef debug
119  LogDebug("BscSim") << "BscSD : Initialize called for " << name << std::endl;
120 #endif
121 
123  if (hcID<0)
124  hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
125  HCE->AddHitsCollection(hcID, theHC);
126 
127  tsID = -2;
128  primID = -2;
129 
131 }
132 
133 
134 bool BscSD::ProcessHits(G4Step * aStep, G4TouchableHistory * ) {
135 
136  if (aStep == NULL) {
137  return true;
138  } else {
139  GetStepInfo(aStep);
140  // LogDebug("BscSim") << edeposit <<std::endl;
141 
142  //AZ
143 #ifdef debug
144  LogDebug("BscSim") << "BscSD : number of hits = " << theHC->entries() << std::endl;
145 #endif
146 
147  if (HitExists() == false && edeposit>0. ){
148  CreateNewHit();
149  return true;
150  }
151  }
152  return true;
153 }
154 
155 void BscSD::GetStepInfo(G4Step* aStep) {
156 
157  preStepPoint = aStep->GetPreStepPoint();
158  postStepPoint= aStep->GetPostStepPoint();
159  theTrack = aStep->GetTrack();
160  hitPoint = preStepPoint->GetPosition();
161  currentPV = preStepPoint->GetPhysicalVolume();
162  hitPointExit = postStepPoint->GetPosition();
163 
164  hitPointLocal = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
165  hitPointLocalExit = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPointExit);
166 
167 
168  G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
169  LogDebug("BscSim") << " BscSD :particleType = " << theTrack->GetDefinition()->GetParticleName() <<std::endl;
170  if (particleCode == emPDG ||
171  particleCode == epPDG ||
172  particleCode == gammaPDG ) {
173  edepositEM = getEnergyDeposit(aStep); edepositHAD = 0.;
174  } else {
175  edepositEM = 0.; edepositHAD = getEnergyDeposit(aStep);
176  }
177  edeposit = aStep->GetTotalEnergyDeposit();
178  tSlice = (postStepPoint->GetGlobalTime() )/nanosecond;
179  tSliceID = (int) tSlice;
180  unitID = setDetUnitId(aStep);
181 #ifdef debug
182  LogDebug("BscSim") << "unitID=" << unitID <<std::endl;
183 #endif
184  primaryID = theTrack->GetTrackID();
185  // Position = hitPoint;
186  Pabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
187  Tof = aStep->GetPostStepPoint()->GetGlobalTime()/nanosecond;
188  Eloss = aStep->GetTotalEnergyDeposit()/GeV;
189  ParticleType = theTrack->GetDefinition()->GetPDGEncoding();
190  ThetaAtEntry = aStep->GetPreStepPoint()->GetPosition().theta()/deg;
191  PhiAtEntry = aStep->GetPreStepPoint()->GetPosition().phi()/deg;
192 
193  ParentId = theTrack->GetParentID();
194  Vx = theTrack->GetVertexPosition().x();
195  Vy = theTrack->GetVertexPosition().y();
196  Vz = theTrack->GetVertexPosition().z();
197  X = hitPoint.x();
198  Y = hitPoint.y();
199  Z = hitPoint.z();
200 }
201 
202 uint32_t BscSD::setDetUnitId(G4Step * aStep) {
203 
204  return (numberingScheme == 0 ? 0 : numberingScheme->getUnitID(aStep));
205 }
206 
207 
209  if (primaryID<1) {
210  edm::LogWarning("BscSim") << "***** BscSD error: primaryID = "
211  << primaryID
212  << " maybe detector name changed";
213  }
214 
215  // Update if in the same detector, time-slice and for same track
216  // if (primaryID == primID && tSliceID == tsID && unitID==previousUnitID) {
217  if (tSliceID == tsID && unitID==previousUnitID) {
218  //AZ:
219  UpdateHit();
220  return true;
221  }
222  // Reset entry point for new primary
223  if (primaryID != primID)
225 
226  //look in the HitContainer whether a hit with the same primID, unitID,
227  //tSliceID already exists:
228 
229  G4bool found = false;
230 
231  // LogDebug("BscSim") << "BscSD: HCollection= " << theHC->entries() <<std::endl;
232 
233  for (int j=0; j<theHC->entries()&&!found; j++) {
234  BscG4Hit* aPreviousHit = (*theHC)[j];
235  if (aPreviousHit->getTrackID() == primaryID &&
236  aPreviousHit->getTimeSliceID() == tSliceID &&
237  aPreviousHit->getUnitID() == unitID ) {
238  //AZ:
239  currentHit = aPreviousHit;
240  found = true;
241  }
242  }
243 
244  if (found) {
245  //AZ:
246  UpdateHit();
247  return true;
248  } else {
249  return false;
250  }
251 }
252 
253 
255 
258  incidentEnergy = preStepPoint->GetKineticEnergy();
259 
260 }
261 
262 
264 
265  if (primID<0) return;
266  if (hit == 0 ) {
267  edm::LogWarning("BscSim") << "BscSD: hit to be stored is NULL !!";
268  return;
269  }
270 
271  theHC->insert( hit );
272 }
273 
274 
276 
277 #ifdef debug
278  LogDebug("BscSim") << "BscSD CreateNewHit for"
279  << " PV " << currentPV->GetName()
280  << " PVid = " << currentPV->GetCopyNo()
281  << " Unit " << unitID <<std::endl;
282  LogDebug("BscSim") << " primary " << primaryID
283  << " time slice " << tSliceID
284  << " For Track " << theTrack->GetTrackID()
285  << " which is a " << theTrack->GetDefinition()->GetParticleName();
286 
287  if (theTrack->GetTrackID()==1) {
288  LogDebug("BscSim") << " of energy " << theTrack->GetTotalEnergy();
289  } else {
290  LogDebug("BscSim") << " daughter of part. " << theTrack->GetParentID();
291  }
292 
293  LogDebug("BscSim") << " and created by " ;
294  if (theTrack->GetCreatorProcess()!=NULL)
295  LogDebug("BscSim") << theTrack->GetCreatorProcess()->GetProcessName() ;
296  else
297  LogDebug("BscSim") << "NO process";
298  LogDebug("BscSim") << std::endl;
299 #endif
300 
301  currentHit = new BscG4Hit;
306 
313 
315 
318 
320  currentHit->setVx(Vx);
321  currentHit->setVy(Vy);
322  currentHit->setVz(Vz);
323 
324  currentHit->setX(X);
325  currentHit->setY(Y);
326  currentHit->setZ(Z);
327 
328  UpdateHit();
329 
331 }
332 
333 
335 
336  if (Eloss > 0.) {
338 
339 #ifdef debug
340  LogDebug("BscSim") << "updateHit: add eloss " << Eloss <<std::endl;
341  LogDebug("BscSim") << "CurrentHit=" << currentHit
342  << ", PostStepPoint=" << postStepPoint->GetPosition();
343 #endif
344  //AZ
346  }
347 
348  // buffer for next steps:
349  tsID = tSliceID;
350  primID = primaryID;
352 }
353 
354 
355 G4ThreeVector BscSD::SetToLocal(const G4ThreeVector& global){
356 
357  const G4VTouchable* touch= preStepPoint->GetTouchable();
358  theEntryPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
359  return theEntryPoint;
360 }
361 
362 
363 G4ThreeVector BscSD::SetToLocalExit(const G4ThreeVector& globalPoint){
364 
365  const G4VTouchable* touch= postStepPoint->GetTouchable();
366  theExitPoint = touch->GetHistory()->GetTopTransform().TransformPoint(globalPoint);
367  return theExitPoint;
368 }
369 
370 
371 void BscSD::EndOfEvent(G4HCofThisEvent* ) {
372 
373  // here we loop over transient hits and make them persistent
374  for (int j=0; j<theHC->entries(); j++) {
375  //AZ:
376  BscG4Hit* aHit = (*theHC)[j];
377  LogDebug("BscSim") << "hit number" << j << "unit ID = "<<aHit->getUnitID()<< "\n";
378  LogDebug("BscSim") << "entry z " << aHit->getEntry().z()<< "\n";
379  LogDebug("BscSim") << "entr theta " << aHit->getThetaAtEntry()<< "\n";
380 
381  Local3DPoint locExitPoint(0,0,0);
382  Local3DPoint locEntryPoint(aHit->getEntry().x(),
383  aHit->getEntry().y(),
384  aHit->getEntry().z());
385  slave->processHits(PSimHit(locEntryPoint,locExitPoint,
386  aHit->getPabs(),
387  aHit->getTof(),
388  aHit->getEnergyLoss(),
389  aHit->getParticleType(),
390  aHit->getUnitID(),
391  aHit->getTrackID(),
392  aHit->getThetaAtEntry(),
393  aHit->getPhiAtEntry()));
394  }
395  Summarize();
396 }
397 
398 
400 }
401 
402 
403 void BscSD::clear() {
404 }
405 
406 
408 }
409 
410 
412  LogDebug("BscSim") << "BscSD: Collection " << theHC->GetName() << "\n";
413  theHC->PrintAllHits();
414 }
415 
416 
418  if (slave->name() == n) c=slave->hits();
419 }
420 
421 void BscSD::update (const BeginOfEvent * i) {
422  LogDebug("BscSim") << " Dispatched BeginOfEvent for " << GetName()
423  << " !" ;
424  clearHits();
425  eventno = (*i)()->GetEventID();
426 }
427 
428 void BscSD::update(const BeginOfRun *) {
429 
430  G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
431  G4String particleName;
432  emPDG = theParticleTable->FindParticle(particleName="e-")->GetPDGEncoding();
433  epPDG = theParticleTable->FindParticle(particleName="e+")->GetPDGEncoding();
434  gammaPDG = theParticleTable->FindParticle(particleName="gamma")->GetPDGEncoding();
435 
436 }
437 
438 void BscSD::update (const ::EndOfEvent*) {
439 }
440 
442  //AZ:
443  slave->Initialize();
444 }
445 
446 std::vector<std::string> BscSD::getNames(){
447  std::vector<std::string> temp;
448  temp.push_back(slave->name());
449  return temp;
450 }
451 
#define LogDebug(id)
void setTof(float e)
Definition: BSCG4Hit.cc:154
T getParameter(std::string const &) const
G4int gammaPDG
Definition: BscSD.h:167
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
G4int emPDG
Definition: BscSD.h:165
void fillHits(edm::PSimHitContainer &, std::string use)
Definition: BscSD.cc:417
virtual ~BscSD()
Definition: BscSD.cc:104
std::vector< std::string > logicalNames(std::string &readoutName)
G4int epPDG
Definition: BscSD.h:166
void setEntry(const G4ThreeVector &xyz)
Definition: BSCG4Hit.cc:116
const double GeV
Definition: MathUtil.h:16
void Summarize()
Definition: BscSD.cc:399
uint32_t previousUnitID
Definition: BscSD.h:131
float getPhiAtEntry() const
Definition: BSCG4Hit.cc:159
void setVx(float p)
Definition: BSCG4Hit.cc:177
float Y
Definition: BscSD.h:154
TrackingSlaveSD * slave
Definition: BscSD.h:110
G4ThreeVector SetToLocal(const G4ThreeVector &global)
Definition: BscSD.cc:355
int getParticleType() const
Definition: BSCG4Hit.cc:151
float PhiAtEntry
Definition: BscSD.h:150
G4int getTrackID() const
Definition: BSCG4Hit.cc:133
G4bool HitExists()
Definition: BscSD.cc:208
float Vx
Definition: BscSD.h:153
float edeposit
Definition: BscSD.h:137
G4int tSliceID
Definition: BscSD.h:132
float Tof
Definition: BscSD.h:145
std::string name() const
virtual double getEnergyDeposit(G4Step *step)
Definition: BscSD.cc:113
void setTrackID(int i)
Definition: BSCG4Hit.cc:134
void setEnergyLoss(float e)
Definition: BSCG4Hit.cc:155
int ParentId
Definition: BscSD.h:152
void setParentId(int p)
Definition: BSCG4Hit.cc:174
void GetStepInfo(G4Step *aStep)
Definition: BscSD.cc:155
float X
Definition: BscSD.h:154
void CreateNewHit()
Definition: BscSD.cc:275
virtual bool ProcessHits(G4Step *, G4TouchableHistory *)
Definition: BscSD.cc:134
float edepositEM
Definition: BscSD.h:164
#define NULL
Definition: scimark2.h:8
void setUnitID(unsigned int i)
Definition: BSCG4Hit.cc:137
float Eloss
Definition: BscSD.h:146
float getTof() const
Definition: BSCG4Hit.cc:149
float getPabs() const
Definition: BSCG4Hit.cc:148
G4int primID
Definition: BscSD.h:118
BscG4HitCollection * theHC
Definition: BscSD.h:123
void setVy(float p)
Definition: BSCG4Hit.cc:180
void setY(float t)
Definition: BSCG4Hit.cc:168
G4ThreeVector hitPointLocalExit
Definition: BscSD.h:143
type of data representation of DDCompactView
Definition: DDCompactView.h:77
unsigned int getUnitID() const
Definition: BSCG4Hit.cc:136
virtual void Initialize(G4HCofThisEvent *HCE)
Definition: BscSD.cc:117
G4int tsID
Definition: BscSD.h:126
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
virtual unsigned int getUnitID(const G4Step *aStep) const
virtual uint32_t setDetUnitId(G4Step *)
Definition: BscSD.cc:202
void UpdateHit()
Definition: BscSD.cc:334
int getTimeSliceID() const
Definition: BSCG4Hit.cc:141
G4VPhysicalVolume * currentPV
Definition: BscSD.h:129
void setThetaAtEntry(float t)
Definition: BSCG4Hit.cc:161
void setX(float t)
Definition: BSCG4Hit.cc:165
void setPabs(float e)
Definition: BSCG4Hit.cc:153
void setEntryLocalP(const G4ThreeVector &xyz1)
Definition: BSCG4Hit.cc:119
std::vector< PSimHit > & hits()
uint32_t unitID
Definition: BscSD.h:131
G4int primaryID
Definition: BscSD.h:132
std::string const collectionName[nCollections]
Definition: Collections.h:45
virtual void Initialize()
void setPhiAtEntry(float f)
Definition: BSCG4Hit.cc:162
G4StepPoint * preStepPoint
Definition: BscSD.h:135
virtual void clearHits()
Definition: BscSD.cc:441
ParticleType
std::string name
Definition: BscSD.h:121
float ThetaAtEntry
Definition: BscSD.h:149
int j
Definition: DBlmapReader.cc:9
G4int hcID
Definition: BscSD.h:122
float Z
Definition: BscSD.h:154
float edepositHAD
Definition: BscSD.h:164
virtual void DrawAll()
Definition: BscSD.cc:407
std::vector< std::string > getNames()
Definition: BscSD.cc:446
G4ThreeVector hitPoint
Definition: BscSD.h:139
float getThetaAtEntry() const
Definition: BSCG4Hit.cc:158
BscSD(std::string, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: BscSD.cc:48
virtual void clear()
Definition: BscSD.cc:403
void addEnergyDeposit(double em, double hd)
Definition: BSCG4Hit.cc:143
void setTimeSlice(double d)
Definition: BSCG4Hit.cc:140
void setParticleType(short i)
Definition: BSCG4Hit.cc:156
void StoreHit(BscG4Hit *)
Definition: BscSD.cc:263
void setVz(float p)
Definition: BSCG4Hit.cc:183
G4ThreeVector theEntryPoint
Definition: BscSD.h:114
G4StepPoint * postStepPoint
Definition: BscSD.h:136
G4ThreeVector hitPointExit
Definition: BscSD.h:141
BscNumberingScheme * numberingScheme
Definition: BscSD.h:111
float incidentEnergy
Definition: BscSD.h:117
G4ThreeVector exitPoint
Definition: BscSD.h:113
G4ThreeVector theExitPoint
Definition: BscSD.h:115
G4ThreeVector SetToLocalExit(const G4ThreeVector &globalPoint)
Definition: BscSD.cc:363
G4Track * theTrack
Definition: BscSD.h:128
virtual void AssignSD(std::string &vname)
int eventno
Definition: BscSD.h:160
G4double tSlice
Definition: BscSD.h:133
G4ThreeVector getEntry() const
Definition: BSCG4Hit.cc:115
void setZ(float t)
Definition: BSCG4Hit.cc:171
virtual void PrintAll()
Definition: BscSD.cc:411
G4ThreeVector entrancePoint
Definition: BscSD.h:113
void update(const BeginOfRun *)
This routine will be called when the appropriate signal arrives.
Definition: BscSD.cc:428
virtual bool processHits(const PSimHit &)
std::vector< PSimHit > PSimHitContainer
void setExitLocalP(const G4ThreeVector &xyz1)
Definition: BSCG4Hit.cc:122
float Vz
Definition: BscSD.h:153
G4THitsCollection< BscG4Hit > BscG4HitCollection
void setIncidentEnergy(double e)
Definition: BSCG4Hit.cc:131
virtual void EndOfEvent(G4HCofThisEvent *eventHC)
Definition: BscSD.cc:371
float Pabs
Definition: BscSD.h:144
float Vy
Definition: BscSD.h:153
float getEnergyLoss() const
Definition: BSCG4Hit.cc:150
G4ThreeVector hitPointLocal
Definition: BscSD.h:142
BscG4Hit * currentHit
Definition: BscSD.h:127
void ResetForNewPrimary()
Definition: BscSD.cc:254