CMS 3D CMS Logo

PltSD.cc
Go to the documentation of this file.
2 
7 
8 #include "G4Step.hh"
9 #include "G4StepPoint.hh"
10 #include "G4Track.hh"
11 #include "G4ThreeVector.hh"
12 
13 #include "CLHEP/Units/GlobalSystemOfUnits.h"
14 #include "CLHEP/Units/GlobalPhysicalConstants.h"
15 
16 #include <iostream>
17 
19  const SensitiveDetectorCatalog & clg,
20  edm::ParameterSet const & p, const SimTrackManager* manager)
21  : TimingSD(name, cpv, clg, p, manager) {
22 
23  edm::ParameterSet m_TrackerSD = p.getParameter<edm::ParameterSet>("PltSD");
24  energyCut = m_TrackerSD.getParameter<double>("EnergyThresholdForPersistencyInGeV")*CLHEP::GeV; //default must be 0.5
25  energyHistoryCut = m_TrackerSD.getParameter<double>("EnergyThresholdForHistoryInGeV")*CLHEP::GeV;//default must be 0.05
26 
27  setCuts(energyCut, energyHistoryCut);
28 }
29 
31 }
32 
33 uint32_t PltSD::setDetUnitId(const G4Step * aStep ) {
34 
35  unsigned int detId = 0;
36 
37  LogDebug("PltSD")<< " DetID = "<<detId;
38 
39  //Find number of levels
40  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
41  int level = 0;
42  if (touch) level = ((touch->GetHistoryDepth())+1);
43 
44  //Get name and copy numbers
45  if ( level > 1 ) {
46  //some debugging with the names
47  G4String sensorName = touch->GetVolume(2)->GetName();
48  G4String telName = touch->GetVolume(3)->GetName();
49  G4String volumeName = touch->GetVolume(4)->GetName();
50  if ( sensorName != "PLTSensorPlane" )
51  std::cout << " PltSD::setDetUnitId -w- Sensor name not PLTSensorPlane " << std::endl;
52  if ( telName != "Telescope" )
53  std::cout << " PltSD::setDetUnitId -w- Telescope name not Telescope " << std::endl;
54  if ( volumeName != "PLT" )
55  std::cout << " PltSD::setDetUnitId -w- Volume name not PLT " << std::endl;
56 
57  //Get the information about which telescope, plane, row/column was hit
58  int columnNum = touch->GetReplicaNumber(0);
59  int rowNum = touch->GetReplicaNumber(1);
60  int sensorNum = touch->GetReplicaNumber(2);
61  int telNum = touch->GetReplicaNumber(3);
62  //temp stores the PLTBCM volume the hit occured in (i.e. was the hit on the + or -z side?)
63  int temp = touch->GetReplicaNumber(5);
64  //convert to the PLT hit id standard
65  int pltNum;
66  if (temp == 2) pltNum = 0;
67  else pltNum = 1;
68 
69  //correct the telescope numbers on the -z side to have the same naming convention in phi as the +z side
70  if (pltNum == 0){
71  if (telNum == 0){
72  telNum = 7;
73  }
74  else if (telNum == 1){
75  telNum = 6;
76  }
77  else if (telNum == 2){
78  telNum = 5;
79  }
80  else if (telNum == 3){
81  telNum = 4;
82  }
83  else if (telNum == 4){
84  telNum = 3;
85  }
86  else if (telNum == 5){
87  telNum = 2;
88  }
89  else if (telNum == 6){
90  telNum = 1;
91  }
92  else if (telNum == 7){
93  telNum = 0;
94  }
95  }
96  //the PLT is divided into sets of telescopes on the + and -x sides
97  int halfCarriageNum = -1;
98 
99  //If the telescope is on the -x side of the carriage, halfCarriageNum=0. If on the +x side, it is = 1.
100  if(telNum == 0 || telNum == 1 || telNum == 2 || telNum == 3)
101  halfCarriageNum = 0;
102  else
103  halfCarriageNum = 1;
104  //correct the telescope numbers of the +x half-carriage to range from 0 to 3
105  if(halfCarriageNum == 1){
106  if(telNum == 4){
107  telNum = 0;
108  }
109  else if (telNum == 5){
110  telNum = 1;
111  }
112  else if (telNum == 6){
113  telNum = 2;
114  }
115  else if (telNum == 7){
116  telNum = 3;
117  }
118  }
119  //Define unique detId for each pixel. See https://twiki.cern.ch/twiki/bin/viewauth/CMS/PLTSimulationGuide for more information
120  detId = 10000000*pltNum+1000000*halfCarriageNum+100000*telNum+10000*sensorNum+100*rowNum+columnNum;
121  //std::cout << "Hit Recorded at " << "plt:" << pltNum << " hc:" << halfCarriageNum << " tel:" << telNum << " plane:" << sensorNum << std::endl;
122  }
123  return detId;
124 }
125 
126 bool PltSD::checkHit(const G4Step*, BscG4Hit* hit) {
127 
128  // 50 micron are allowed between the exit
129  // point of the current hit and the entry point of the new hit
130  static const float tolerance2 = (float)(0.0025 * CLHEP::mm * CLHEP::mm);
131  return ((hit->getExitLocalP()-getLocalEntryPoint()).mag2()<tolerance2);
132 }
#define LogDebug(id)
T getParameter(std::string const &) const
const double GeV
Definition: MathUtil.h:16
PltSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: PltSD.cc:18
const G4ThreeVector & getLocalEntryPoint() const
Definition: TimingSD.h:69
void setCuts(double eCut, double historyCut)
Definition: TimingSD.cc:87
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:83
bool checkHit(const G4Step *, BscG4Hit *) override
Definition: PltSD.cc:126
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
const G4ThreeVector & getExitLocalP() const
Definition: BscG4Hit.h:38
double energyHistoryCut
Definition: PltSD.h:29
~PltSD() override
Definition: PltSD.cc:30
double energyCut
Definition: PltSD.h:28
uint32_t setDetUnitId(const G4Step *) override
Definition: PltSD.cc:33