CMS 3D CMS Logo

MtdSD.cc
Go to the documentation of this file.
2 
5 
9 
10 #include "G4Track.hh"
11 #include "G4Step.hh"
12 #include "G4StepPoint.hh"
13 
14 #include <iostream>
15 
16 //#define EDM_ML_DEBUG
17 //-------------------------------------------------------------------
19  const SensitiveDetectorCatalog& clg,
20  edm::ParameterSet const& p,
21  const SimTrackManager* manager)
22  : TimingSD(name, clg, manager), numberingScheme(nullptr) {
23  //Parameters
24  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("MtdSD");
25  int verbn = m_p.getUntrackedParameter<int>("Verbosity");
26 
27  SetVerboseLevel(verbn);
28 
29  MTDNumberingScheme* scheme = nullptr;
30  if (name == "FastTimerHitsBarrel") {
31  scheme = dynamic_cast<MTDNumberingScheme*>(new BTLNumberingScheme());
32  isBTL = true;
33  } else if (name == "FastTimerHitsEndcap") {
34  scheme = dynamic_cast<MTDNumberingScheme*>(new ETLNumberingScheme());
35  isETL = true;
36  } else {
37  scheme = nullptr;
38  edm::LogWarning("MtdSim") << "MtdSD: ReadoutName not supported";
39  }
40  if (scheme)
42 
43  double newTimeFactor = 1. / m_p.getParameter<double>("TimeSliceUnit");
44  edm::LogVerbatim("MtdSim") << "New time factor = " << newTimeFactor;
45  setTimeFactor(newTimeFactor);
46 
47  edm::LogVerbatim("MtdSim") << "MtdSD: Instantiation completed for " << name;
48 }
49 
51 
52 uint32_t MtdSD::setDetUnitId(const G4Step* aStep) {
53  if (numberingScheme == nullptr) {
54  return MTDDetId();
55  } else {
56  getBaseNumber(aStep);
57 #ifdef EDM_ML_DEBUG
58  edm::LogVerbatim("MtdSim") << "DetId = " << numberingScheme->getUnitID(theBaseNumber);
59 #endif
61  }
62 }
63 
65  if (scheme != nullptr) {
66  edm::LogVerbatim("MtdSim") << "MtdSD: updates numbering scheme for " << GetName();
67  if (numberingScheme)
68  delete numberingScheme;
70  }
71 }
72 
73 void MtdSD::getBaseNumber(const G4Step* aStep) {
75  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
76  int theSize = touch->GetHistoryDepth() + 1;
77  if (theBaseNumber.getCapacity() < theSize)
78  theBaseNumber.setSize(theSize);
79  //Get name and copy numbers
80  if (theSize > 1) {
81 #ifdef EDM_ML_DEBUG
82  edm::LogVerbatim("MtdSim") << "Building MTD basenumber:";
83 #endif
84  for (int ii = 0; ii < theSize; ii++) {
85  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(), touch->GetReplicaNumber(ii));
86 #ifdef EDM_ML_DEBUG
87  edm::LogVerbatim("MtdSim") << "MtdSD::getBaseNumber(): Adding level " << ii << ": "
88  << touch->GetVolume(ii)->GetName() << "[" << touch->GetReplicaNumber(ii) << "]";
89 #endif
90  }
91  }
92 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void setTimeFactor(double)
Definition: TimingSD.cc:79
bool isBTL
Definition: MtdSD.h:33
bool isETL
Definition: MtdSD.h:34
void addLevel(const std::string_view name, const int copyNumber)
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
uint32_t setDetUnitId(const G4Step *) override
Definition: MtdSD.cc:52
T getUntrackedParameter(std::string const &, T const &) const
MTDBaseNumber theBaseNumber
Definition: MtdSD.h:32
void getBaseNumber(const G4Step *)
Definition: MtdSD.cc:73
void setNumberingScheme(MTDNumberingScheme *)
Definition: MtdSD.cc:64
ii
Definition: cuy.py:589
MtdSD(const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: MtdSD.cc:18
MTDNumberingScheme * numberingScheme
Definition: MtdSD.h:31
virtual uint32_t getUnitID(const MTDBaseNumber &baseNumber) const =0
static TrackerG4SimHitNumberingScheme & numberingScheme(const GeometricDet &det)
Log< level::Warning, false > LogWarning
void setSize(const int size)
~MtdSD() override
Definition: MtdSD.cc:50