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("MtdSim") << "MtdSD: ReadoutName not supported";
58  }
59  if (scheme) setNumberingScheme(scheme);
60 
61  double newTimeFactor = 1./m_p.getParameter<double>("TimeSliceUnit");
62  edm::LogInfo("MtdSim") << "New time factor = " << newTimeFactor;
63  setTimeFactor(newTimeFactor);
64 
65  edm::LogVerbatim("MtdSim") << "MtdSD: Instantiation completed for "
66  << name << " of type " << type;
67 }
68 
70 }
71 
72 uint32_t MtdSD::setDetUnitId(const G4Step * aStep) {
73  if (numberingScheme == nullptr) {
74  return MTDDetId();
75  } else {
76  getBaseNumber(aStep);
78  }
79 }
80 
81 std::vector<double> MtdSD::getDDDArray(const std::string & str,
82  const DDsvalues_type & sv) {
83 
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
90  << " bins " << nval << " < 1 ==> illegal";
91  throw cms::Exception("DDException") << "MtdSD: cannot get array " << str;
92  }
93  return fvec;
94  } else {
95  edm::LogError("MtdSim") << "MtdSD: cannot get array " << str;
96  throw cms::Exception("DDException") << "MtdSD: cannot get array " << str;
97  }
98 }
99 
101  if (scheme != nullptr) {
102  edm::LogInfo("MtdSim") << "MtdSD: updates numbering scheme for "
103  << GetName();
104  if (numberingScheme) delete numberingScheme;
105  numberingScheme = scheme;
106  }
107 }
108 
109 void MtdSD::getBaseNumber(const G4Step* aStep) {
110 
112  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
113  int theSize = touch->GetHistoryDepth()+1;
114  if ( theBaseNumber.getCapacity() < theSize ) theBaseNumber.setSize(theSize);
115  //Get name and copy numbers
116  if ( theSize > 1 ) {
117  for (int ii = 0; ii < theSize ; ii++) {
118  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(),touch->GetReplicaNumber(ii));
119 #ifdef EDM_ML_DEBUG
120  edm::LogInfo("MtdSim") << "MtdSD::getBaseNumber(): Adding level " << ii
121  << ": " << touch->GetVolume(ii)->GetName() << "["
122  << touch->GetReplicaNumber(ii) << "]";
123 #endif
124  }
125  }
126 }
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:80
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:72
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
void getBaseNumber(const G4Step *)
Definition: MtdSD.cc:109
void setNumberingScheme(MTDNumberingScheme *)
Definition: MtdSD.cc:100
ii
Definition: cuy.py:590
MTDNumberingScheme * numberingScheme
Definition: MtdSD.h:37
void setSize(const int &size)
DDsvalues_type mergedSpecifics() const
void addLevel(const std::string &name, const int &copyNumber)
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
Definition: DDsvalues.h:12
bool firstChild()
set the current node to the first child ...
#define str(s)
~MtdSD() override
Definition: MtdSD.cc:69
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: MtdSD.cc:81