CMS 3D CMS Logo

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