CMS 3D CMS Logo

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