CMS 3D CMS Logo

MtdSD.cc
Go to the documentation of this file.
1 //#define EDM_ML_DEBUG
2 
4 
7 
11 
12 #include "G4Track.hh"
13 #include "G4Step.hh"
14 #include "G4StepPoint.hh"
15 
16 #include <iostream>
17 
18 //-------------------------------------------------------------------
20  const SensitiveDetectorCatalog& clg,
21  edm::ParameterSet const& p,
22  const SimTrackManager* manager)
23  : TimingSD(name, clg, manager), numberingScheme(nullptr) {
24  //Parameters
25  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("MtdSD");
26  int verbn = m_p.getUntrackedParameter<int>("Verbosity");
27 
28  SetVerboseLevel(verbn);
29 
30  MTDNumberingScheme* scheme = nullptr;
31  if (name == "FastTimerHitsBarrel") {
32  scheme = dynamic_cast<MTDNumberingScheme*>(new BTLNumberingScheme());
33  isBTL = true;
34  } else if (name == "FastTimerHitsEndcap") {
35  scheme = dynamic_cast<MTDNumberingScheme*>(new ETLNumberingScheme());
36  isETL = true;
37  } else {
38  scheme = nullptr;
39  edm::LogWarning("MtdSim") << "MtdSD: ReadoutName not supported";
40  }
41  if (scheme)
43 
44  energyCut = m_p.getParameter<double>("EnergyThresholdForPersistencyInGeV") * CLHEP::GeV; //default must be 0.5
45  energyHistoryCut = m_p.getParameter<double>("EnergyThresholdForHistoryInGeV") * CLHEP::GeV; //default must be 0.05
46 
48 
49  double newTimeFactor = 1. / m_p.getParameter<double>("TimeSliceUnit");
50  edm::LogVerbatim("MtdSim") << "New time factor = " << newTimeFactor;
51  setTimeFactor(newTimeFactor);
52 
53  edm::LogVerbatim("MtdSim") << "MtdSD: Instantiation completed for " << name;
54 }
55 
57 
58 uint32_t MtdSD::setDetUnitId(const G4Step* aStep) {
59  if (numberingScheme == nullptr) {
60  return MTDDetId();
61  } else {
62  getBaseNumber(aStep);
63 #ifdef EDM_ML_DEBUG
64  edm::LogVerbatim("MtdSim") << "DetId = " << numberingScheme->getUnitID(theBaseNumber);
65 #endif
67  }
68 }
69 
71  if (scheme != nullptr) {
72  edm::LogVerbatim("MtdSim") << "MtdSD: updates numbering scheme for " << GetName();
73  if (numberingScheme)
74  delete numberingScheme;
76  }
77 }
78 
79 void MtdSD::getBaseNumber(const G4Step* aStep) {
81  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
82  int theSize = touch->GetHistoryDepth() + 1;
83  if (theBaseNumber.getCapacity() < theSize)
84  theBaseNumber.setSize(theSize);
85  //Get name and copy numbers
86  if (theSize > 1) {
87 #ifdef EDM_ML_DEBUG
88  edm::LogVerbatim("MtdSim") << "Building MTD basenumber:";
89 #endif
90  for (int ii = 0; ii < theSize; ii++) {
91  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(), touch->GetReplicaNumber(ii));
92 #ifdef EDM_ML_DEBUG
93  edm::LogVerbatim("MtdSim") << "MtdSD::getBaseNumber(): Adding level " << ii << ": "
94  << touch->GetVolume(ii)->GetName() << "[" << touch->GetReplicaNumber(ii) << "]";
95 #endif
96  }
97  }
98 }
99 
100 int MtdSD::getTrackID(const G4Track* aTrack) {
101  int theID = aTrack->GetTrackID();
102  TrackInformation* trkInfo = cmsTrackInformation(aTrack);
103  const G4String& rname = aTrack->GetVolume()->GetLogicalVolume()->GetRegion()->GetName();
104  if (trkInfo != nullptr) {
105 #ifdef EDM_ML_DEBUG
106  trkInfo->Print();
107 #endif
108  if (rname == "FastTimerRegionSensBTL") {
109  if (trkInfo->isBTLdaughter()) {
110  theID = trkInfo->idAtBTLentrance();
111  }
112 #ifdef EDM_ML_DEBUG
113  edm::LogVerbatim("MtdSim") << "BTL Track ID: " << trkInfo->idAtBTLentrance() << ":" << theID;
114 #endif
115  } else if (rname == "FastTimerRegionSensETL") {
116  theID = trkInfo->getIDonCaloSurface();
117 #ifdef EDM_ML_DEBUG
118  edm::LogVerbatim("MtdSim") << "ETL Track ID: " << trkInfo->getIDonCaloSurface() << ":" << theID;
119 #endif
120  } else {
121  throw cms::Exception("MtdSDError") << "MtdSD called in incorrect region " << rname;
122  }
123  } else {
124 #ifdef EDM_ML_DEBUG
125  edm::LogWarning("MtdSim") << "MtdSD: Problem with primaryID **** set by force to TkID **** " << theID;
126 #endif
127  }
128  return theID;
129 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void setTimeFactor(double)
Definition: TimingSD.cc:79
bool isBTL
Definition: MtdSD.h:39
bool isETL
Definition: MtdSD.h:40
bool isBTLdaughter() const
void addLevel(const std::string_view name, const int copyNumber)
void setCuts(double eCut, double historyCut)
Definition: TimingSD.cc:89
int getTrackID(const G4Track *) override
Definition: MtdSD.cc:100
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
uint32_t setDetUnitId(const G4Step *) override
Definition: MtdSD.cc:58
T getUntrackedParameter(std::string const &, T const &) const
double energyCut
Definition: MtdSD.h:31
MTDBaseNumber theBaseNumber
Definition: MtdSD.h:38
double energyHistoryCut
Definition: MtdSD.h:32
void getBaseNumber(const G4Step *)
Definition: MtdSD.cc:79
void setNumberingScheme(MTDNumberingScheme *)
Definition: MtdSD.cc:70
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
int getIDonCaloSurface() const
ii
Definition: cuy.py:589
MtdSD(const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: MtdSD.cc:19
void Print() const override
MTDNumberingScheme * numberingScheme
Definition: MtdSD.h:37
const G4String rname[NREG]
virtual uint32_t getUnitID(const MTDBaseNumber &baseNumber) const =0
static TrackerG4SimHitNumberingScheme & numberingScheme(const GeometricDet &det)
int idAtBTLentrance() const
Log< level::Warning, false > LogWarning
void setSize(const int size)
~MtdSD() override
Definition: MtdSD.cc:56