CMS 3D CMS Logo

Bcm1fSD.cc
Go to the documentation of this file.
2 
5 
9 
12 
15 
18 
19 #include "G4Track.hh"
20 #include "G4SDManager.hh"
21 #include "G4VProcess.hh"
22 #include "G4EventManager.hh"
23 #include "G4Event.hh"
24 #include "G4VProcess.hh"
25 
26 #include <string>
27 #include <vector>
28 #include <iostream>
29 
30 #include <CLHEP/Units/SystemOfUnits.h>
31 
33  const SensitiveDetectorCatalog& clg,
34  edm::ParameterSet const& p,
35  const SimTrackManager* manager)
36  : TimingSD(name, clg, manager) {
37  edm::ParameterSet m_TrackerSD = p.getParameter<edm::ParameterSet>("Bcm1fSD");
38  energyCut =
39  m_TrackerSD.getParameter<double>("EnergyThresholdForPersistencyInGeV") * CLHEP::GeV; //default must be 0.5 (?)
41  m_TrackerSD.getParameter<double>("EnergyThresholdForHistoryInGeV") * CLHEP::GeV; //default must be 0.05 (?)
42 
44 }
45 
47 
48 uint32_t Bcm1fSD::setDetUnitId(const G4Step* aStep) {
49  uint32_t detId = 0;
50 
51  //Find number of levels
52  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
53  int level = (touch) ? ((touch->GetHistoryDepth()) + 1) : 0;
54 
55  //Get name and copy numbers
56  if (level > 1) {
57  std::string sensorName = ForwardName::getName(touch->GetVolume(0)->GetName());
58  std::string diamondName = ForwardName::getName(touch->GetVolume(1)->GetName());
59  std::string detectorName = ForwardName::getName(touch->GetVolume(2)->GetName());
60  std::string volumeName = ForwardName::getName(touch->GetVolume(3)->GetName());
61 
62  if (sensorName != "BCM1FSensor") {
63  edm::LogWarning("ForwardSim") << "Bcm1fSD::setDetUnitId -w- Sensor name " << sensorName << " not BCM1FSensor ";
64  }
65  if (detectorName != "BCM1F") {
66  edm::LogWarning("ForwardSim") << " Bcm1fSD::setDetUnitId -w- Detector name " << detectorName << " not BCM1F ";
67  }
68  int sensorNo = touch->GetReplicaNumber(0);
69  int diamondNo = touch->GetReplicaNumber(1);
70  int volumeNo = touch->GetReplicaNumber(3);
71 
72  // Detector ID definition
73  // detId = XYYZ
74  // X = volume, 1: +Z, 2: -Z
75  // YY = diamond, 01-12, 12: phi = 90 deg, numbering clockwise when looking from the IP
76  // Z = sensor, 1 or 2, clockwise when looking from the IP
77 
78  detId = 1000 * volumeNo + 10 * diamondNo + sensorNo;
79  }
80  return detId;
81 }
82 
83 bool Bcm1fSD::checkHit(const G4Step*, BscG4Hit* hit) {
84  // 50 micron are allowed between the exit
85  // point of the current hit and the entry point of the new hit
86  static const float tolerance2 = (float)(0.0025 * CLHEP::mm * CLHEP::mm);
87  return ((hit->getExitLocalP() - getLocalEntryPoint()).mag2() < tolerance2);
88 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
~Bcm1fSD() override
Definition: Bcm1fSD.cc:46
void setCuts(double eCut, double historyCut)
Definition: TimingSD.cc:89
float energyHistoryCut
Definition: Bcm1fSD.h:25
uint32_t setDetUnitId(const G4Step *) override
Definition: Bcm1fSD.cc:48
Bcm1fSD(const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: Bcm1fSD.cc:32
const G4ThreeVector & getLocalEntryPoint() const
Definition: TimingSD.h:62
bool checkHit(const G4Step *, BscG4Hit *) override
Definition: Bcm1fSD.cc:83
float energyCut
Definition: Bcm1fSD.h:24
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
std::string getName(const G4String &)
Definition: ForwardName.cc:3
Log< level::Warning, false > LogWarning