CMS 3D CMS Logo

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