CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PLTSensitiveDetector.cc
Go to the documentation of this file.
2 
5 
9 
12 
14 
19 
20 #include "G4Track.hh"
21 #include "G4SDManager.hh"
22 #include "G4VProcess.hh"
23 #include "G4EventManager.hh"
24 #include "G4Event.hh"
25 #include "G4VProcess.hh"
26 
27 #include <string>
28 #include <vector>
29 #include <iostream>
30 
31 
33  const DDCompactView & cpv,
35  edm::ParameterSet const & p,
36  const SimTrackManager* manager) :
37  SensitiveTkDetector(name, cpv, clg, p), myName(name), mySimHit(0),
38  oldVolume(0), lastId(0), lastTrack(0), eventno(0) {
39 
40  edm::ParameterSet m_TrackerSD = p.getParameter<edm::ParameterSet>("PltSD");
41  energyCut = m_TrackerSD.getParameter<double>("EnergyThresholdForPersistencyInGeV")*GeV; //default must be 0.5
42  energyHistoryCut = m_TrackerSD.getParameter<double>("EnergyThresholdForHistoryInGeV")*GeV;//default must be 0.05
43 
44  edm::LogInfo("PLTSensitiveDetector") <<"Criteria for Saving Tracker SimTracks: \n "
45  <<" History: "<<energyHistoryCut<< " MeV ; Persistency: "<< energyCut<<" MeV\n"
46  <<" Constructing a PLTSensitiveDetector with ";
47 
48  slave = new TrackingSlaveSD(name);
49 
50  // Now attach the right detectors (LogicalVolumes) to me
51  std::vector<std::string> lvNames = clg.logicalNames(name);
52  this->Register();
53  for (std::vector<std::string>::iterator it = lvNames.begin();
54  it != lvNames.end(); it++) {
55  edm::LogInfo("PLTSensitiveDetector")<< name << " attaching LV " << *it;
56  this->AssignSD(*it);
57  }
58 
61 }
62 
64  delete slave;
66  delete myG4TrackToParticleID;
67 }
68 
69 bool PLTSensitiveDetector::ProcessHits(G4Step * aStep, G4TouchableHistory *) {
70 
71  LogDebug("PLTSensitiveDetector") << " Entering a new Step "
72  << aStep->GetTotalEnergyDeposit() << " "
73  << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
74 
75  G4Track * gTrack = aStep->GetTrack();
76  if ((unsigned int)(gTrack->GetTrackID()) != lastTrack) {
77 
78  if (gTrack->GetKineticEnergy() > energyCut){
80  info->storeTrack(true);
81  }
82  if (gTrack->GetKineticEnergy() > energyHistoryCut){
84  info->putInHistory();
85  }
86  }
87 
88  if (aStep->GetTotalEnergyDeposit()>0.) {
89  if (newHit(aStep) == true) {
90  sendHit();
91  createHit(aStep);
92  } else {
93  updateHit(aStep);
94  }
95  return true;
96  }
97  return false;
98 }
99 
101 
102  unsigned int detId = 0;
103 
104  LogDebug("PLTSensitiveDetector")<< " DetID = "<<detId;
105 
106  return detId;
107 }
108 
109 void PLTSensitiveDetector::EndOfEvent(G4HCofThisEvent *) {
110 
111  LogDebug("PLTSensitiveDetector")<< " Saving the last hit in a ROU " << myName;
112 
113  if (mySimHit == 0) return;
114  sendHit();
115 }
116 
118  if (slave->name() == n) c=slave->hits();
119 }
120 
122  if (mySimHit == 0) return;
123  LogDebug("PLTSensitiveDetector") << " Storing PSimHit: " << pname << " " << mySimHit->detUnitId()
124  << " " << mySimHit->trackId() << " " << mySimHit->energyLoss()
125  << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
126 
128 
129  // clean up
130  delete mySimHit;
131  mySimHit = 0;
132  lastTrack = 0;
133  lastId = 0;
134 }
135 
136 void PLTSensitiveDetector::updateHit(G4Step * aStep) {
137 
139  float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
140  mySimHit->setExitPoint(theExitPoint);
141  LogDebug("PLTSensitiveDetector")<< " Before : " << mySimHit->energyLoss();
142  mySimHit->addEnergyLoss(theEnergyLoss);
144 
145  LogDebug("PLTSensitiveDetector") << " Updating: new exitpoint " << pname << " "
146  << theExitPoint << " new energy loss " << theEnergyLoss
147  << "\n Updated PSimHit: " << mySimHit->detUnitId()
148  << " " << mySimHit->trackId()
149  << " " << mySimHit->energyLoss() << " "
150  << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
151 }
152 
153 bool PLTSensitiveDetector::newHit(G4Step * aStep) {
154 
155  G4Track * theTrack = aStep->GetTrack();
156  uint32_t theDetUnitId = setDetUnitId(aStep);
157  unsigned int theTrackID = theTrack->GetTrackID();
158 
159  LogDebug("PLTSensitiveDetector") << " OLD (d,t) = (" << lastId << "," << lastTrack
160  << "), new = (" << theDetUnitId << "," << theTrackID << ") return "
161  << ((theTrackID == lastTrack) && (lastId == theDetUnitId));
162  if ((mySimHit != 0) && (theTrackID == lastTrack) && (lastId == theDetUnitId) && closeHit(aStep))
163  return false;
164  return true;
165 }
166 
167 bool PLTSensitiveDetector::closeHit(G4Step * aStep) {
168 
169  if (mySimHit == 0) return false;
170  const float tolerance = 0.05 * mm; // 50 micron are allowed between the exit
171  // point of the current hit and the entry point of the new hit
173  LogDebug("PLTSensitiveDetector")<< " closeHit: distance = " << (mySimHit->exitPoint()-theEntryPoint).mag();
174 
175  if ((mySimHit->exitPoint()-theEntryPoint).mag()<tolerance) return true;
176  return false;
177 }
178 
179 void PLTSensitiveDetector::createHit(G4Step * aStep) {
180 
181  if (mySimHit != 0) {
182  delete mySimHit;
183  mySimHit=0;
184  }
185 
186  G4Track * theTrack = aStep->GetTrack();
187  G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume();
188 
191 
192  float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
193  float theTof = aStep->GetPreStepPoint()->GetGlobalTime()/nanosecond;
194  float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
195  int theParticleType = myG4TrackToParticleID->particleID(theTrack);
196  uint32_t theDetUnitId = setDetUnitId(aStep);
197 
200  pname = theTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
201 
202  unsigned int theTrackID = theTrack->GetTrackID();
203 
204  G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection();
205  // convert it to local frame
206  G4ThreeVector lmd = ((G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory()->GetTopTransform().TransformAxis(gmd);
208  float theThetaAtEntry = lnmd.theta();
209  float thePhiAtEntry = lnmd.phi();
210 
211  mySimHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof,
212  theEnergyLoss,theParticleType,theDetUnitId,
213  theTrackID,theThetaAtEntry,thePhiAtEntry,
214  theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess()));
215 
216  LogDebug("PLTSensitiveDetector") << " Created PSimHit: " << pname << " "
217  << mySimHit->detUnitId() << " " << mySimHit->trackId()
218  << " " << mySimHit->energyLoss() << " "
219  << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
220  lastId = theDetUnitId;
221  lastTrack = theTrackID;
222  oldVolume = v;
223 }
224 
226 
228 
229  clearHits();
230  eventno = (*i)()->GetEventID();
231  mySimHit = 0;
232 }
233 
235 
236  const G4Track* gTrack = (*bot)();
237  pname = gTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
238 }
239 
241  slave->Initialize();
242 }
243 
245  G4VUserTrackInformation* temp = gTrack->GetUserInformation();
246  if (temp == 0){
247  edm::LogError("PLTSensitiveDetector") <<" ERROR: no G4VUserTrackInformation available";
248  abort();
249  }else{
250  TrackInformation* info = dynamic_cast<TrackInformation*>(temp);
251  if (info == 0){
252  edm::LogError("PLTSensitiveDetector") <<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation";
253  abort();
254  }
255  return info;
256  }
257 }
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
UpdatablePSimHit * mySimHit
virtual void createHit(G4Step *)
std::vector< std::string > logicalNames(std::string &readoutName)
G4TrackToParticleID * myG4TrackToParticleID
bool storeTrack() const
void update(const BeginOfEvent *)
This routine will be called when the appropriate signal arrives.
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:68
TrackInformation * getOrCreateTrackInformation(const G4Track *)
void setExitPoint(const Local3DPoint &exit)
G4ProcessTypeEnumerator * theG4ProcessTypeEnumerator
G4VPhysicalVolume * oldVolume
type of data representation of DDCompactView
Definition: DDCompactView.h:77
virtual bool newHit(G4Step *)
Geom::Theta< T > theta() const
Definition: PV3DBase.h:74
virtual bool closeHit(G4Step *)
std::vector< PSimHit > & hits()
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
virtual void Initialize()
TrackingSlaveSD * slave
virtual void updateHit(G4Step *)
void fillHits(edm::PSimHitContainer &, std::string use)
int particleID(const G4Track *)
virtual uint32_t setDetUnitId(G4Step *)
void addEnergyLoss(float eloss)
unsigned int processId(const G4VProcess *)
virtual void AssignSD(std::string &vname)
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
virtual void EndOfEvent(G4HCofThisEvent *)
unsigned int trackId() const
Definition: PSimHit.h:102
virtual bool ProcessHits(G4Step *, G4TouchableHistory *)
virtual bool processHits(const PSimHit &)
PLTSensitiveDetector(std::string, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
std::vector< PSimHit > PSimHitContainer
Local3DPoint FinalStepPosition(G4Step *s, coordinates)
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
mathSSE::Vec4< T > v
Local3DPoint InitialStepPosition(G4Step *s, coordinates)
unsigned int detUnitId() const
Definition: PSimHit.h:93