CMS 3D CMS Logo

TkAccumulatingSensitiveDetector.cc
Go to the documentation of this file.
3 
8 
11 
18 
22 
25 
27 
28 #include "G4Track.hh"
29 #include "G4SDManager.hh"
30 #include "G4VProcess.hh"
31 #include "G4EventManager.hh"
32 #include "G4Event.hh"
33 
34 #include "G4SystemOfUnits.hh"
35 
36 
37 #include <string>
38 #include <vector>
39 #include <iostream>
40 
42 
43 //#define FAKEFRAMEROTATION
44 //#define DUMPPROCESSES
45 
46 #ifdef DUMPPROCESSES
47 #include "G4VProcess.hh"
48 #endif
49 
50 using std::cout;
51 using std::endl;
52 using std::vector;
53 using std::string;
54 
55 static
58 {
59  static thread_local TrackerG4SimHitNumberingScheme s_scheme(cpv, det);
60  return s_scheme;
61 }
62 
63 
65  const DDCompactView & cpv,
66  const SensitiveDetectorCatalog & clg,
67  edm::ParameterSet const & p,
68  const SimTrackManager* manager) :
69  SensitiveTkDetector(name, cpv, clg, p), myName(name), myRotation(nullptr), mySimHit(nullptr),theManager(manager),
70  oldVolume(nullptr), lastId(0), lastTrack(0), eventno(0) ,rTracker(1200.*mm),zTracker(3000.*mm),
71  numberingScheme_(nullptr)
72 {
73 
74  edm::ParameterSet m_TrackerSD = p.getParameter<edm::ParameterSet>("TrackerSD");
75  allowZeroEnergyLoss = m_TrackerSD.getParameter<bool>("ZeroEnergyLoss");
76  neverAccumulate = m_TrackerSD.getParameter<bool>("NeverAccumulate");
77  printHits = m_TrackerSD.getParameter<bool>("PrintHits");
78  theSigma = m_TrackerSD.getParameter<double>("ElectronicSigmaInNanoSeconds");
79  energyCut = m_TrackerSD.getParameter<double>("EnergyThresholdForPersistencyInGeV")*GeV; //default must be 0.5
80  energyHistoryCut = m_TrackerSD.getParameter<double>("EnergyThresholdForHistoryInGeV")*GeV;//default must be 0.05
81 
82  edm::LogInfo("TrackerSimInfo") <<"Criteria for Saving Tracker SimTracks: ";
83  edm::LogInfo("TrackerSimInfo")<<" History: "<<energyHistoryCut<< " MeV ; Persistency: "<< energyCut<<" MeV ";
84  edm::LogInfo("TrackerSimInfo")<<" Constructing a TkAccumulatingSensitiveDetector with ";
85 
86 #ifndef FAKEFRAMEROTATION
87  // No Rotation given in input, automagically choose one based upon the name
88  if (name.find("TrackerHits") != string::npos)
89  {
90  edm::LogInfo("TrackerSimInfo")<<" TkAccumulatingSensitiveDetector: using TrackerFrameRotation for "<<myName;
92  }
93  // Just in case (test beam etc)
94  if (myRotation == nullptr)
95  {
96  edm::LogInfo("TrackerSimInfo")<<" TkAccumulatingSensitiveDetector: using StandardFrameRotation for "<<myName;
98  }
99 #else
101  edm::LogWarning("TrackerSimInfo")<< " WARNING - Using FakeFrameRotation in TkAccumulatingSensitiveDetector;";
102 #endif
103 
104  slaveLowTof = new TrackingSlaveSD(name+"LowTof");
105  slaveHighTof = new TrackingSlaveSD(name+"HighTof");
106 
107  // Now attach the right detectors (LogicalVolumes) to me
108  const vector<string>& lvNames = clg.logicalNames(name);
109  this->Register();
110  for (vector<string>::const_iterator it = lvNames.begin(); it != lvNames.end(); it++)
111  {
112  edm::LogInfo("TrackerSimInfo")<< name << " attaching LV " << *it;
113  this->AssignSD(*it);
114  }
115 
118 }
119 
121 {
122  delete slaveLowTof;
123  delete slaveHighTof;
125  delete myG4TrackToParticleID;
126 }
127 
128 
129 bool TkAccumulatingSensitiveDetector::ProcessHits(G4Step * aStep, G4TouchableHistory * ROhist)
130 {
131 
132  LogDebug("TrackerSimDebug")<< " Entering a new Step " << aStep->GetTotalEnergyDeposit() << " "
133  << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
134 
135  //TimeMe t1(theTimer, false);
136  if (aStep->GetTotalEnergyDeposit()>0. || allowZeroEnergyLoss == true)
137  {
138  if (newHit(aStep) == true)
139  {
140  sendHit();
141  createHit(aStep);
142  }
143  else
144  {
145  updateHit(aStep);
146  }
147  return true;
148  }
149  return false;
150 }
151 
153 {
154  unsigned int detId = numberingScheme_->g4ToNumberingScheme(s->GetPreStepPoint()->GetTouchable());
155 
156  LogDebug("TrackerSimDebug")<< " DetID = "<<detId;
157 
158  return detId;
159 }
160 
161 
163 {
164 
165  float threshold = 3. * theSigma;
166  if (tof < threshold) return 1;
167  return 2;
168 }
169 
171 {
172  if (myRotation !=nullptr) return myRotation->transformPoint(in,v);
173  return (in);
174 }
175 
176 
178  const G4Track* gTrack = (*bot)();
179 #ifdef DUMPPROCESSES
180  edm::LogInfo("TrackerSimInfo") <<" -> process creator pointer "<<gTrack->GetCreatorProcess();
181  if (gTrack->GetCreatorProcess())
182  edm::LogInfo("TrackerSimInfo")<<" -> PROCESS CREATOR : "<<gTrack->GetCreatorProcess()->GetProcessName();
183 #endif
184 
185 
186  //
187  //Position
188  //
189  const G4ThreeVector pos = gTrack->GetPosition();
190  //LogDebug("TrackerSimDebug")<<" ENERGY MeV "<<gTrack->GetKineticEnergy()<<" Energy Cut" << energyCut;
191  //LogDebug("TrackerSimDebug")<<" TOTAL ENERGY "<<gTrack->GetTotalEnergy();
192  //LogDebug("TrackerSimDebug")<<" WEIGHT "<<gTrack->GetWeight();
193 
194  //
195  // Check if in Tracker Volume
196  //
197  if (pos.perp() < rTracker &&std::fabs(pos.z()) < zTracker){
198  //
199  // inside the Tracker
200  //
201  LogDebug("TrackerSimDebug")<<" INSIDE TRACKER";
202 
203  if (gTrack->GetKineticEnergy() > energyCut){
205  //LogDebug("TrackerSimDebug")<<" POINTER "<<info;
206  //LogDebug("TrackerSimDebug")<<" track inside the tracker selected for STORE";
207  //LogDebug("TrackerSimDebug")<<"Track ID (persistent track) = "<<gTrack->GetTrackID();
208 
209  info->storeTrack(true);
210  }
211  //
212  // Save History?
213  //
214  if (gTrack->GetKineticEnergy() > energyHistoryCut){
216  info->putInHistory();
217  //LogDebug("TrackerSimDebug")<<" POINTER "<<info;
218  //LogDebug("TrackerSimDebug")<<" track inside the tracker selected for HISTORY";
219  //LogDebug("TrackerSimDebug")<<"Track ID (history track) = "<<gTrack->GetTrackID();
220  }
221 
222  }
223 }
224 
225 
226 
228 {
229  if (mySimHit == nullptr) return;
230  LogDebug("TrackerSimDebug")<< " Storing PSimHit: " << pname << " " << mySimHit->detUnitId()
231  << " " << mySimHit->trackId() << " " << mySimHit->energyLoss()
232  << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
233  if (printHits)
234  {
235  TkSimHitPrinter thePrinter("TkHitPositionOSCAR.dat");
236  thePrinter.startNewSimHit(myName,oldVolume->GetLogicalVolume()->GetName(),
238  thePrinter.printLocal(mySimHit->entryPoint(),mySimHit->exitPoint());
241  thePrinter.printMomentumOfTrack(mySimHit->pabs(),pname,
243  thePrinter.printGlobalMomentum(px,py,pz);
244  }
245 
246  if (tofBin(mySimHit->timeOfFlight()) == 1)
247  slaveLowTof->processHits(*mySimHit); // implicit conversion (slicing) to PSimHit!!!
248  else
249  slaveHighTof->processHits(*mySimHit); // implicit conversion (slicing) to PSimHit!!!
250  //
251  // clean up
252  delete mySimHit;
253  mySimHit = nullptr;
254  lastTrack = 0;
255  lastId = 0;
256 }
257 
259 {
260  if (mySimHit != nullptr)
261  {
262  delete mySimHit;
263  mySimHit=nullptr;
264  }
265 
266  G4Track * theTrack = aStep->GetTrack();
267 
268  G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume();
269  Local3DPoint theEntryPoint;
270  Local3DPoint theExitPoint =
272 
273  //
274  // Check particle type - for gamma and neutral hadrons energy deposition
275  // should be local (V.I.)
276  //
277  if(0.0 == theTrack->GetDefinition()->GetPDGCharge()) {
278  theEntryPoint = theExitPoint;
279  } else {
281  }
282 
283  //
284  // This allows to send he skipEvent if it is outside!
285  //
286  checkExitPoint(theExitPoint);
287  float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
288  float theTof = aStep->GetPreStepPoint()->GetGlobalTime()/nanosecond;
289  float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
290  int theParticleType = myG4TrackToParticleID->particleID(theTrack);
291  uint32_t theDetUnitId = setDetUnitId(aStep);
292 
293  //
294  // Check on particle charge is not applied because these points are not stored
295  // in hits (V.I.)
296  //
299 
300  pname = theTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
301 
302  if (theDetUnitId == 0)
303  {
304  edm::LogError("TrackerSimInfo") << " Error: theDetUnitId is not valid.";
305  abort();
306  }
307  unsigned int theTrackID = theTrack->GetTrackID();
308 
309  // To whom assign the Hit?
310  // First iteration: if the track is to be stored, use the current number;
311  // otherwise, get to the mother
312  unsigned int theTrackIDInsideTheSimHit=theTrackID;
313 
314 
315  G4VUserTrackInformation * info = theTrack->GetUserInformation();
316  if (info == nullptr) edm::LogError("TrackerSimInfo")<< " Error: no UserInformation available ";
317  else
318  {
319  TrackInformation * temp = dynamic_cast<TrackInformation* >(info);
320  if (temp ==nullptr) edm::LogError("TrackerSimInfo")<< " Error:G4VUserTrackInformation is not a TrackInformation.";
321  else {
322  if (temp->storeTrack() == false) {
323  // Go to the mother!
324  LogDebug("TrackerSimDebug")<< " TkAccumulatingSensitiveDetector:createHit(): setting the TrackID from "
325  << theTrackIDInsideTheSimHit;
326  theTrackIDInsideTheSimHit = theTrack->GetParentID();
327  LogDebug("TrackerSimDebug")<< " to the mother one " << theTrackIDInsideTheSimHit << " " << theEnergyLoss;
328  }
329  else {
330  LogDebug("TrackerSimDebug")<< " TkAccumulatingSensitiveDetector:createHit(): leaving the current TrackID "
331  << theTrackIDInsideTheSimHit;
332  }
333  }
334  }
335 
336  px = aStep->GetPreStepPoint()->GetMomentum().x()/GeV;
337  py = aStep->GetPreStepPoint()->GetMomentum().y()/GeV;
338  pz = aStep->GetPreStepPoint()->GetMomentum().z()/GeV;
339 
340  G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection();
341  // convert it to local frame
342  G4ThreeVector lmd = ((G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory()
343  ->GetTopTransform().TransformAxis(gmd);
345  float theThetaAtEntry = lnmd.theta();
346  float thePhiAtEntry = lnmd.phi();
347 
348  mySimHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof,
349  theEnergyLoss,theParticleType,theDetUnitId,
350  theTrackIDInsideTheSimHit,theThetaAtEntry,thePhiAtEntry,
351  theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess()));
352  LogDebug("TrackerSimDebug")<< " Created PSimHit: " << pname << " " << mySimHit->detUnitId() << " " << mySimHit->trackId()
353  << " " << mySimHit->energyLoss() << " " << mySimHit->entryPoint()
354  << " " << mySimHit->exitPoint();
355  lastId = theDetUnitId;
356  lastTrack = theTrackID;
357  oldVolume = v;
358 }
359 
361 {
362  G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume();
364  //
365  // This allows to send he skipEvent if it is outside!
366  //
367  checkExitPoint(theExitPoint);
368  float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
369  mySimHit->setExitPoint(theExitPoint);
370  LogDebug("TrackerSimDebug")<< " Before : " << mySimHit->energyLoss();
371  mySimHit->addEnergyLoss(theEnergyLoss);
373 
374  LogDebug("TrackerSimDebug")<< " Updating: new exitpoint " << pname << " " << theExitPoint
375  << " new energy loss " << theEnergyLoss;
376  LogDebug("TrackerSimDebug")<< " Updated PSimHit: " << mySimHit->detUnitId() << " " << mySimHit->trackId()
377  << " " << mySimHit->energyLoss() << " " << mySimHit->entryPoint()
378  << " " << mySimHit->exitPoint();
379 }
380 
382 {
383  if (neverAccumulate == true) return true;
384  G4Track * theTrack = aStep->GetTrack();
385 
386  // for neutral particles do not merge hits (V.I.)
387  if(0.0 == theTrack->GetDefinition()->GetPDGCharge()) return true;
388 
389  uint32_t theDetUnitId = setDetUnitId(aStep);
390  unsigned int theTrackID = theTrack->GetTrackID();
391 
392  LogDebug("TrackerSimDebug")<< " OLD (d,t) = (" << lastId << "," << lastTrack
393  << "), new = (" << theDetUnitId << "," << theTrackID << ") return "
394  << ((theTrackID == lastTrack) && (lastId == theDetUnitId));
395  if ((mySimHit != nullptr) && (theTrackID == lastTrack) && (lastId == theDetUnitId) && closeHit(aStep))
396  return false;
397  return true;
398 }
399 
401 {
402  if (mySimHit == nullptr) return false;
403  const float tolerance = 0.05 * mm; // 50 micron are allowed between the exit
404  // point of the current hit and the entry point of the new hit
405  G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume();
407  LogDebug("TrackerSimDebug")<< " closeHit: distance = " << (mySimHit->exitPoint()-theEntryPoint).mag();
408 
409  if ((mySimHit->exitPoint()-theEntryPoint).mag()<tolerance) return true;
410  return false;
411 }
412 
414 {
415 
416  LogDebug("TrackerSimDebug")<< " Saving the last hit in a ROU " << myName;
417 
418  if (mySimHit == nullptr) return;
419  sendHit();
420 }
421 
423 {
424  clearHits();
425  eventno = (*i)()->GetEventID();
426  mySimHit = nullptr;
427 }
428 
430 {
432  const edm::EventSetup* es = (*i)();
433  es->get<IdealGeometryRecord>().get( pDD );
434 
436  es->get<IdealGeometryRecord>().get(pView);
437 
438  numberingScheme_=&(numberingScheme(*pView,*pDD));
439 }
440 
442 {
445 }
446 
448 {
449  double z = p.z();
450  if (std::abs(z)<0.3*mm) return;
451  bool sendExc = false;
452  //static SimpleConfigurable<bool> sendExc(false,"TrackerSim:ThrowOnBadHits");
453  edm::LogWarning("TrackerSimInfo")<< " ************ Hit outside the detector ; Local z " << z
454  << "; skipping event = " << sendExc;
455  if (sendExc == true)
456  {
457  G4EventManager::GetEventManager()->AbortCurrentEvent();
458  G4EventManager::GetEventManager()->GetNonconstCurrentEvent()->SetEventAborted();
459  }
460 }
461 
463  G4VUserTrackInformation* temp = gTrack->GetUserInformation();
464  if (temp == nullptr){
465  edm::LogError("TrackerSimInfo") <<" ERROR: no G4VUserTrackInformation available";
466  abort();
467  }else{
468  TrackInformation* info = dynamic_cast<TrackInformation*>(temp);
469  if (info ==nullptr){
470  edm::LogError("TrackerSimInfo")<<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation";
471  abort();
472  }
473  return info;
474  }
475 }
476 
478  //
479  // do it once for low, once for High
480  //
481 
482  if (slaveLowTof->name() == n) c=slaveLowTof->hits();
483  if (slaveHighTof->name() == n) c=slaveHighTof->hits();
484 
485 }
486 
488  std::vector<string> temp;
489  temp.push_back(slaveLowTof->name());
490  temp.push_back(slaveHighTof->name());
491  return temp;
492 }
#define LogDebug(id)
T getParameter(std::string const &) const
TrackInformation * getOrCreateTrackInformation(const G4Track *)
static const TGPicture * info(bool iBackgroundIsBlack)
const double GeV
Definition: MathUtil.h:16
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
bool storeTrack() const
void fillHits(edm::PSimHitContainer &, std::string use) override
Local3DPoint toOrcaRef(Local3DPoint, G4VPhysicalVolume *)
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::string name() const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
void printMomentumOfTrack(float, std::string, int sign) const
void setExitPoint(const Local3DPoint &exit)
bool ProcessHits(G4Step *, G4TouchableHistory *) override
void printHitData(float, float) const
const std::vector< std::string > & logicalNames(const std::string &readoutName) const
#define nullptr
type of data representation of DDCompactView
Definition: DDCompactView.h:90
void printGlobalMomentum(float, float, float) const
virtual void AssignSD(const std::string &vname)
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
unsigned int processId(const G4VProcess *p)
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
std::vector< PSimHit > & hits()
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point)
float timeOfFlight() const
Definition: PSimHit.h:69
virtual void Initialize()
unsigned int g4ToNumberingScheme(const G4VTouchable *)
T z() const
Definition: PV3DBase.h:64
void printGlobal(Local3DPoint, Local3DPoint) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:63
static int particleID(const G4Track *)
void addEnergyLoss(float eloss)
const T & get() const
Definition: EventSetup.h:55
TrackerG4SimHitNumberingScheme * numberingScheme_
std::vector< std::string > getNames() override
void startNewSimHit(std::string, std::string, int, int, int)
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
void printLocal(Local3DPoint, Local3DPoint) const
unsigned int trackId() const
Definition: PSimHit.h:102
virtual Local3DPoint transformPoint(Local3DPoint &, G4VPhysicalVolume *) const =0
virtual bool processHits(const PSimHit &)
std::vector< PSimHit > PSimHitContainer
void EndOfEvent(G4HCofThisEvent *) override
Local3DPoint FinalStepPosition(G4Step *s, coordinates)
int getPropagationSign(Local3DPoint, Local3DPoint)
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
Local3DPoint InitialStepPosition(G4Step *s, coordinates)
TkAccumulatingSensitiveDetector(std::string, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
unsigned int detUnitId() const
Definition: PSimHit.h:93
G4ProcessTypeEnumerator * theG4ProcessTypeEnumerator