CMS 3D CMS Logo

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