CMS 3D CMS Logo

MtdSD.cc
Go to the documentation of this file.
2 
9 
13 
15 
19 
20 #include "G4Track.hh"
21 #include "G4Step.hh"
22 #include "G4StepPoint.hh"
23 
24 #include <iostream>
25 
26 //#define EDM_ML_DEBUG
27 //-------------------------------------------------------------------
29  const edm::EventSetup& es,
30  const SensitiveDetectorCatalog& clg,
31  edm::ParameterSet const& p,
32  const SimTrackManager* manager)
33  : TimingSD(name, es, clg, p, manager), numberingScheme(nullptr) {
34  //Parameters
36  int verbn = m_p.getUntrackedParameter<int>("Verbosity");
37 
38  SetVerboseLevel(verbn);
39 
41  es.get<IdealGeometryRecord>().get(cpv);
42 
43  std::string attribute = "ReadOutName";
44  DDSpecificsMatchesValueFilter filter{DDValue(attribute, name, 0)};
45  DDFilteredView fv(*cpv, filter);
46  fv.firstChild();
48  std::vector<int> temp = dbl_to_int(getDDDArray("Type", sv));
49  int type = temp[0];
50 
51  MTDNumberingScheme* scheme = nullptr;
52  if (name == "FastTimerHitsBarrel") {
53  scheme = dynamic_cast<MTDNumberingScheme*>(new BTLNumberingScheme());
54  isBTL = true;
55  } else if (name == "FastTimerHitsEndcap") {
56  scheme = dynamic_cast<MTDNumberingScheme*>(new ETLNumberingScheme());
57  isETL = true;
58  } else {
59  scheme = nullptr;
60  edm::LogWarning("MtdSim") << "MtdSD: ReadoutName not supported";
61  }
62  if (scheme)
63  setNumberingScheme(scheme);
64 
65  double newTimeFactor = 1. / m_p.getParameter<double>("TimeSliceUnit");
66  edm::LogInfo("MtdSim") << "New time factor = " << newTimeFactor;
67  setTimeFactor(newTimeFactor);
68 
69  edm::LogVerbatim("MtdSim") << "MtdSD: Instantiation completed for " << name << " of type " << type;
70 }
71 
73 
74 uint32_t MtdSD::setDetUnitId(const G4Step* aStep) {
75  if (numberingScheme == nullptr) {
76  return MTDDetId();
77  } else {
78  getBaseNumber(aStep);
80  }
81 }
82 
83 std::vector<double> MtdSD::getDDDArray(const std::string& str, const DDsvalues_type& sv) {
84  DDValue value(str);
85  if (DDfetch(&sv, value)) {
86  const std::vector<double>& fvec = value.doubles();
87  int nval = fvec.size();
88  if (nval < 1) {
89  edm::LogError("MtdSim") << "MtdSD : # of " << str << " bins " << nval << " < 1 ==> illegal";
90  throw cms::Exception("DDException") << "MtdSD: cannot get array " << str;
91  }
92  return fvec;
93  } else {
94  edm::LogError("MtdSim") << "MtdSD: cannot get array " << str;
95  throw cms::Exception("DDException") << "MtdSD: cannot get array " << str;
96  }
97 }
98 
100  if (scheme != nullptr) {
101  edm::LogInfo("MtdSim") << "MtdSD: updates numbering scheme for " << GetName();
102  if (numberingScheme)
103  delete numberingScheme;
105  }
106 }
107 
108 void MtdSD::getBaseNumber(const G4Step* aStep) {
110  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
111  int theSize = touch->GetHistoryDepth() + 1;
112  if (theBaseNumber.getCapacity() < theSize)
113  theBaseNumber.setSize(theSize);
114  //Get name and copy numbers
115  if (theSize > 1) {
116  for (int ii = 0; ii < theSize; ii++) {
117  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(), touch->GetReplicaNumber(ii));
118 #ifdef EDM_ML_DEBUG
119  edm::LogInfo("MtdSim") << "MtdSD::getBaseNumber(): Adding level " << ii << ": " << touch->GetVolume(ii)->GetName()
120  << "[" << touch->GetReplicaNumber(ii) << "]";
121 #endif
122  }
123  }
124 }
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const std::vector< double > & doubles() const
a reference to the double-valued values stored in the given instance of DDValue
Definition: DDValue.cc:111
void setTimeFactor(double)
Definition: TimingSD.cc:86
bool isBTL
Definition: MtdSD.h:37
bool isETL
Definition: MtdSD.h:38
#define nullptr
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:79
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
uint32_t setDetUnitId(const G4Step *) override
Definition: MtdSD.cc:74
virtual uint32_t getUnitID(const MTDBaseNumber &baseNumber) const =0
std::vector< int > dbl_to_int(const std::vector< double > &vecdbl)
Converts a std::vector of doubles to a std::vector of int.
Definition: DDutils.h:7
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
Definition: DDsvalues.h:12
MTDBaseNumber theBaseNumber
Definition: MtdSD.h:36
void getBaseNumber(const G4Step *)
Definition: MtdSD.cc:108
void setNumberingScheme(MTDNumberingScheme *)
Definition: MtdSD.cc:99
ii
Definition: cuy.py:590
MTDNumberingScheme * numberingScheme
Definition: MtdSD.h:35
void setSize(const int &size)
DDsvalues_type mergedSpecifics() const
void addLevel(const std::string &name, const int &copyNumber)
MtdSD(const std::string &, const edm::EventSetup &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: MtdSD.cc:28
T get() const
Definition: EventSetup.h:73
bool firstChild()
set the current node to the first child ...
static TrackerG4SimHitNumberingScheme & numberingScheme(const GeometricDet &det)
#define str(s)
~MtdSD() override
Definition: MtdSD.cc:72
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: MtdSD.cc:83