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 SensitiveDetectorCatalog & clg,
30  edm::ParameterSet const & p,
31  const SimTrackManager* manager) :
32  TimingSD(name, cpv, clg, p, manager), numberingScheme(nullptr) {
33 
34  //Parameters
36  int verbn = m_p.getUntrackedParameter<int>("Verbosity");
37 
38  SetVerboseLevel(verbn);
39 
40  std::string attribute = "ReadOutName";
41  DDSpecificsMatchesValueFilter filter{DDValue(attribute,name,0)};
42  DDFilteredView fv(cpv,filter);
43  fv.firstChild();
45  std::vector<int> temp = dbl_to_int(getDDDArray("Type",sv));
46  int type = temp[0];
47 
48  MTDNumberingScheme* scheme=nullptr;
49  if (name == "FastTimerHitsBarrel") {
50  scheme = dynamic_cast<MTDNumberingScheme*>(new BTLNumberingScheme());
51  isBTL=true;
52  } else if (name == "FastTimerHitsEndcap") {
53  scheme = dynamic_cast<MTDNumberingScheme*>(new ETLNumberingScheme());
54  isETL=true;
55  } else {
56  scheme = nullptr;
57  edm::LogWarning("EcalSim") << "MtdSD: ReadoutName not supported";
58  }
59  if (scheme) setNumberingScheme(scheme);
60 
61  setTimeFactor(100.);
62 
63  edm::LogVerbatim("MtdSim") << "MtdSD: Instantiation completed for "
64  << name << " of type " << type;
65 }
66 
68 }
69 
70 uint32_t MtdSD::setDetUnitId(const G4Step * aStep) {
71  if (numberingScheme == nullptr) {
72  return MTDDetId();
73  } else {
74  getBaseNumber(aStep);
76  }
77 }
78 
79 std::vector<double> MtdSD::getDDDArray(const std::string & str,
80  const DDsvalues_type & sv) {
81 
82  DDValue value(str);
83  if (DDfetch(&sv,value)) {
84  const std::vector<double> & fvec = value.doubles();
85  int nval = fvec.size();
86  if (nval < 1) {
87  edm::LogError("MtdSim") << "MtdSD : # of " << str
88  << " bins " << nval << " < 1 ==> illegal";
89  throw cms::Exception("DDException") << "MtdSD: cannot get array " << str;
90  }
91  return fvec;
92  } else {
93  edm::LogError("MtdSim") << "MtdSD: cannot get array " << str;
94  throw cms::Exception("DDException") << "MtdSD: cannot get array " << str;
95  }
96 }
97 
99  if (scheme != nullptr) {
100  edm::LogInfo("MtdSim") << "MtdSD: updates numbering scheme for "
101  << GetName();
102  if (numberingScheme) delete numberingScheme;
103  numberingScheme = scheme;
104  }
105 }
106 
107 void MtdSD::getBaseNumber(const G4Step* aStep) {
108 
110  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
111  int theSize = touch->GetHistoryDepth()+1;
112  if ( theBaseNumber.getCapacity() < theSize ) theBaseNumber.setSize(theSize);
113  //Get name and copy numbers
114  if ( theSize > 1 ) {
115  for (int ii = 0; ii < theSize ; ii++) {
116  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(),touch->GetReplicaNumber(ii));
117 #ifdef EDM_ML_DEBUG
118  edm::LogInfo("MtdSim") << "MtdSD::getBaseNumber(): Adding level " << ii
119  << ": " << 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:140
void setTimeFactor(double)
Definition: TimingSD.cc:77
bool isBTL
Definition: MtdSD.h:39
bool isETL
Definition: MtdSD.h:40
#define nullptr
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:83
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:81
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
uint32_t setDetUnitId(const G4Step *) override
Definition: MtdSD.cc:70
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
MtdSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: MtdSD.cc:28
MTDBaseNumber theBaseNumber
Definition: MtdSD.h:38
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
std::maps an index to a DDValue. The index corresponds to the index assigned to the name of the std::...
Definition: DDsvalues.h:20
void getBaseNumber(const G4Step *)
Definition: MtdSD.cc:107
void setNumberingScheme(MTDNumberingScheme *)
Definition: MtdSD.cc:98
ii
Definition: cuy.py:589
MTDNumberingScheme * numberingScheme
Definition: MtdSD.h:37
void setSize(const int &size)
DDsvalues_type mergedSpecifics() const
void addLevel(const std::string &name, const int &copyNumber)
bool firstChild()
set the current node to the first child ...
#define str(s)
~MtdSD() override
Definition: MtdSD.cc:67
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: MtdSD.cc:79