CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PltSD.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 "CLHEP/Units/GlobalSystemOfUnits.h"
28 #include "CLHEP/Units/GlobalPhysicalConstants.h"
29 
30 #include <string>
31 #include <vector>
32 #include <iostream>
33 
34 
36  const DDCompactView & cpv,
38  edm::ParameterSet const & p,
39  const SimTrackManager* manager) :
40 SensitiveTkDetector(name, cpv, clg, p), myName(name), mySimHit(0),
41 oldVolume(0), lastId(0), lastTrack(0), eventno(0) {
42 
43  edm::ParameterSet m_TrackerSD = p.getParameter<edm::ParameterSet>("PltSD");
44  energyCut = m_TrackerSD.getParameter<double>("EnergyThresholdForPersistencyInGeV")*GeV; //default must be 0.5 (?)
45  energyHistoryCut = m_TrackerSD.getParameter<double>("EnergyThresholdForHistoryInGeV")*GeV;//default must be 0.05 (?)
46 
47  edm::LogInfo("PltSD") <<"Criteria for Saving Tracker SimTracks: \n "
48  <<" History: "<<energyHistoryCut<< " MeV ; Persistency: "<< energyCut<<" MeV\n"
49  <<" Constructing a PltSD with ";
50 
51  slave = new TrackingSlaveSD(name);
52 
53  // Now attach the right detectors (LogicalVolumes) to me
54  std::vector<std::string> lvNames = clg.logicalNames(name);
55  this->Register();
56  for (std::vector<std::string>::iterator it = lvNames.begin();
57  it != lvNames.end(); it++) {
58  edm::LogInfo("PltSD")<< name << " attaching LV " << *it;
59  this->AssignSD(*it);
60  }
61 
64 }
65 
67  delete slave;
69  delete myG4TrackToParticleID;
70 }
71 
72 bool PltSD::ProcessHits(G4Step * aStep, G4TouchableHistory *) {
73 
74  LogDebug("PltSD") << " Entering a new Step "
75  << aStep->GetTotalEnergyDeposit() << " "
76  << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
77 
78  G4Track * gTrack = aStep->GetTrack();
79  if ((unsigned int)(gTrack->GetTrackID()) != lastTrack) {
80 
81  if (gTrack->GetKineticEnergy() > energyCut){
83  info->storeTrack(true);
84  }
85  if (gTrack->GetKineticEnergy() > energyHistoryCut){
87  info->putInHistory();
88  }
89  }
90 
91  if (aStep->GetTotalEnergyDeposit()>0.) {
92  if (newHit(aStep) == true) {
93  sendHit();
94  createHit(aStep);
95  } else {
96  updateHit(aStep);
97  }
98  return true;
99  }
100  return false;
101 }
102 
103 uint32_t PltSD::setDetUnitId(G4Step * aStep ) {
104 
105  unsigned int detId = 0;
106 
107  LogDebug("PltSD")<< " DetID = "<<detId;
108 
109  //Find number of levels
110  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
111  int level = 0;
112  if (touch) level = ((touch->GetHistoryDepth())+1);
113 
114  //Get name and copy numbers
115  if ( level > 1 ) {
116  //some debugging with the names
117  G4String sensorName = touch->GetVolume(2)->GetName();
118  G4String telName = touch->GetVolume(3)->GetName();
119  G4String volumeName = touch->GetVolume(4)->GetName();
120  if ( sensorName != "PLTSensorPlane" )
121  std::cout << " PltSD::setDetUnitId -w- Sensor name not PLTSensorPlane " << std::endl;
122  if ( telName != "Telescope" )
123  std::cout << " PltSD::setDetUnitId -w- Telescope name not Telescope " << std::endl;
124  if ( volumeName != "PLT" )
125  std::cout << " PltSD::setDetUnitId -w- Volume name not PLT " << std::endl;
126 
127  //Get the information about which telescope, plane, row/column was hit
128  int columnNum = touch->GetReplicaNumber(0);
129  int rowNum = touch->GetReplicaNumber(1);
130  int sensorNum = touch->GetReplicaNumber(2);
131  int telNum = touch->GetReplicaNumber(3);
132  //temp stores the PLTBCM volume the hit occured in (i.e. was the hit on the + or -z side?)
133  int temp = touch->GetReplicaNumber(5);
134  //convert to the PLT hit id standard
135  int pltNum;
136  if (temp == 2) pltNum = 0;
137  else pltNum = 1;
138 
139  //correct the telescope numbers on the -z side to have the same naming convention in phi as the +z side
140  if (pltNum == 0){
141  if (telNum == 0){
142  telNum = 7;
143  }
144  else if (telNum == 1){
145  telNum = 6;
146  }
147  else if (telNum == 2){
148  telNum = 5;
149  }
150  else if (telNum == 3){
151  telNum = 4;
152  }
153  else if (telNum == 4){
154  telNum = 3;
155  }
156  else if (telNum == 5){
157  telNum = 2;
158  }
159  else if (telNum == 6){
160  telNum = 1;
161  }
162  else if (telNum == 7){
163  telNum = 0;
164  }
165  }
166  //the PLT is divided into sets of telescopes on the + and -x sides
167  int halfCarriageNum = -1;
168 
169  //If the telescope is on the -x side of the carriage, halfCarriageNum=0. If on the +x side, it is = 1.
170  if(telNum == 0 || telNum == 1 || telNum == 2 || telNum == 3)
171  halfCarriageNum = 0;
172  else
173  halfCarriageNum = 1;
174  //correct the telescope numbers of the +x half-carriage to range from 0 to 3
175  if(halfCarriageNum == 1){
176  if(telNum == 4){
177  telNum = 0;
178  }
179  else if (telNum == 5){
180  telNum = 1;
181  }
182  else if (telNum == 6){
183  telNum = 2;
184  }
185  else if (telNum == 7){
186  telNum = 3;
187  }
188  }
189  //Define unique detId for each pixel. See https://twiki.cern.ch/twiki/bin/viewauth/CMS/PLTSimulationGuide for more information
190  detId = 10000000*pltNum+1000000*halfCarriageNum+100000*telNum+10000*sensorNum+100*rowNum+columnNum;
191  //std::cout << "Hit Recorded at " << "plt:" << pltNum << " hc:" << halfCarriageNum << " tel:" << telNum << " plane:" << sensorNum << std::endl;
192  }
193 
194  return detId;
195 }
196 
197 void PltSD::EndOfEvent(G4HCofThisEvent *) {
198 
199  LogDebug("PltSD")<< " Saving the last hit in a ROU " << myName;
200 
201  if (mySimHit == 0) return;
202  sendHit();
203 }
204 
206  if (slave->name() == n) c=slave->hits();
207 }
208 
210  if (mySimHit == 0) return;
211  LogDebug("PltSD") << " Storing PSimHit: " << pname << " " << mySimHit->detUnitId()
212  << " " << mySimHit->trackId() << " " << mySimHit->energyLoss()
213  << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
214 
216 
217  // clean up
218  delete mySimHit;
219  mySimHit = 0;
220  lastTrack = 0;
221  lastId = 0;
222 }
223 
224 void PltSD::updateHit(G4Step * aStep) {
225 
227  float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
228  mySimHit->setExitPoint(theExitPoint);
229  LogDebug("PltSD")<< " Before : " << mySimHit->energyLoss();
230  mySimHit->addEnergyLoss(theEnergyLoss);
232 
233  LogDebug("PltSD") << " Updating: new exitpoint " << pname << " "
234  << theExitPoint << " new energy loss " << theEnergyLoss
235  << "\n Updated PSimHit: " << mySimHit->detUnitId()
236  << " " << mySimHit->trackId()
237  << " " << mySimHit->energyLoss() << " "
238  << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
239 }
240 
241 bool PltSD::newHit(G4Step * aStep) {
242 
243  G4Track * theTrack = aStep->GetTrack();
244  uint32_t theDetUnitId = setDetUnitId(aStep);
245  unsigned int theTrackID = theTrack->GetTrackID();
246 
247  LogDebug("PltSD") << " OLD (d,t) = (" << lastId << "," << lastTrack
248  << "), new = (" << theDetUnitId << "," << theTrackID << ") return "
249  << ((theTrackID == lastTrack) && (lastId == theDetUnitId));
250  if ((mySimHit != 0) && (theTrackID == lastTrack) && (lastId == theDetUnitId) && closeHit(aStep))
251  return false;
252  return true;
253 }
254 
255 bool PltSD::closeHit(G4Step * aStep) {
256 
257  if (mySimHit == 0) return false;
258  const float tolerance = 0.05 * mm; // 50 micron are allowed between the exit
259  // point of the current hit and the entry point of the new hit
261  LogDebug("PltSD")<< " closeHit: distance = " << (mySimHit->exitPoint()-theEntryPoint).mag();
262 
263  if ((mySimHit->exitPoint()-theEntryPoint).mag()<tolerance) return true;
264  return false;
265 }
266 
267 void PltSD::createHit(G4Step * aStep) {
268 
269  if (mySimHit != 0) {
270  delete mySimHit;
271  mySimHit=0;
272  }
273 
274  G4Track * theTrack = aStep->GetTrack();
275  G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume();
276 
279 
280  float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
281  float theTof = aStep->GetPreStepPoint()->GetGlobalTime()/nanosecond;
282  float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
283  int theParticleType = myG4TrackToParticleID->particleID(theTrack);
284  uint32_t theDetUnitId = setDetUnitId(aStep);
285 
288  pname = theTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
289 
290  unsigned int theTrackID = theTrack->GetTrackID();
291 
292  G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection();
293  // convert it to local frame
294  G4ThreeVector lmd = ((G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory()->GetTopTransform().TransformAxis(gmd);
296  float theThetaAtEntry = lnmd.theta();
297  float thePhiAtEntry = lnmd.phi();
298 
299  mySimHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof,
300  theEnergyLoss,theParticleType,theDetUnitId,
301  theTrackID,theThetaAtEntry,thePhiAtEntry,
302  theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess()));
303 
304  LogDebug("PltSD") << " Created PSimHit: " << pname << " "
305  << mySimHit->detUnitId() << " " << mySimHit->trackId()
306  << " " << mySimHit->energyLoss() << " "
307  << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
308  lastId = theDetUnitId;
309  lastTrack = theTrackID;
310  oldVolume = v;
311 }
312 
313 void PltSD::update(const BeginOfJob * ) { }
314 
315 void PltSD::update(const BeginOfEvent * i) {
316 
317  clearHits();
318  eventno = (*i)()->GetEventID();
319  mySimHit = 0;
320 }
321 
322 void PltSD::update(const BeginOfTrack *bot) {
323 
324  const G4Track* gTrack = (*bot)();
325  pname = gTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
326 }
327 
329  slave->Initialize();
330 }
331 
333  G4VUserTrackInformation* temp = gTrack->GetUserInformation();
334  if (temp == 0){
335  edm::LogError("PltSD") <<" ERROR: no G4VUserTrackInformation available";
336  abort();
337  }else{
338  TrackInformation* info = dynamic_cast<TrackInformation*>(temp);
339  if (info == 0){
340  edm::LogError("PltSD") <<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation";
341  abort();
342  }
343  return info;
344  }
345 }
#define LogDebug(id)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< std::string > logicalNames(std::string &readoutName)
static const TGPicture * info(bool iBackgroundIsBlack)
UpdatablePSimHit * mySimHit
Definition: PltSD.h:67
virtual bool newHit(G4Step *)
Definition: PltSD.cc:241
bool storeTrack() const
std::string myName
Definition: PltSD.h:66
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::string name() const
G4TrackToParticleID * myG4TrackToParticleID
Definition: PltSD.h:65
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
void setExitPoint(const Local3DPoint &exit)
virtual void updateHit(G4Step *)
Definition: PltSD.cc:224
int eventno
Definition: PltSD.h:76
TrackingSlaveSD * slave
Definition: PltSD.h:63
type of data representation of DDCompactView
Definition: DDCompactView.h:77
unsigned int processId(const G4VProcess *p)
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
std::vector< PSimHit > & hits()
float energyCut
Definition: PltSD.h:68
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point)
virtual void Initialize()
virtual void sendHit()
Definition: PltSD.cc:209
float energyHistoryCut
Definition: PltSD.h:69
static int particleID(const G4Track *)
Local3DPoint globalExitPoint
Definition: PltSD.h:72
virtual uint32_t setDetUnitId(G4Step *)
Definition: PltSD.cc:103
void addEnergyLoss(float eloss)
void update(const BeginOfEvent *)
This routine will be called when the appropriate signal arrives.
Definition: PltSD.cc:315
virtual void EndOfEvent(G4HCofThisEvent *)
Definition: PltSD.cc:197
virtual void createHit(G4Step *)
Definition: PltSD.cc:267
virtual void clearHits()
Definition: PltSD.cc:328
virtual void AssignSD(std::string &vname)
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
G4VPhysicalVolume * oldVolume
Definition: PltSD.h:73
virtual bool closeHit(G4Step *)
Definition: PltSD.cc:255
unsigned int trackId() const
Definition: PSimHit.h:102
unsigned int lastTrack
Definition: PltSD.h:75
virtual ~PltSD()
Definition: PltSD.cc:66
virtual bool processHits(const PSimHit &)
tuple cout
Definition: gather_cfg.py:121
std::vector< PSimHit > PSimHitContainer
virtual bool ProcessHits(G4Step *, G4TouchableHistory *)
Definition: PltSD.cc:72
void fillHits(edm::PSimHitContainer &, std::string use)
Definition: PltSD.cc:205
tuple level
Definition: testEve_cfg.py:34
std::string pname
Definition: PltSD.h:77
TrackInformation * getOrCreateTrackInformation(const G4Track *)
Definition: PltSD.cc:332
uint32_t lastId
Definition: PltSD.h:74
Local3DPoint FinalStepPosition(G4Step *s, coordinates)
G4ProcessTypeEnumerator * theG4ProcessTypeEnumerator
Definition: PltSD.h:64
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
Local3DPoint globalEntryPoint
Definition: PltSD.h:71
PltSD(std::string, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: PltSD.cc:35