CMS 3D CMS Logo

FastTimerSD.cc
Go to the documentation of this file.
2 
4 
11 
15 
18 
19 #include <vector>
20 
21 #include "G4Track.hh"
22 #include "G4Step.hh"
23 #include "G4StepPoint.hh"
24 #include "G4VPhysicalVolume.hh"
25 
26 #include <iostream>
27 
28 //#define EDM_ML_DEBUG
29 //-------------------------------------------------------------------
31  const SensitiveDetectorCatalog & clg,
32  edm::ParameterSet const & p,
33  const SimTrackManager* manager) :
34  TimingSD(name, cpv, clg, p, manager), ftcons(nullptr) {
35 
36  //Parameters
37  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("FastTimerSD");
38  int verbn = m_p.getUntrackedParameter<int>("Verbosity");
39 
40  SetVerboseLevel(verbn);
41 
42  std::string attribute = "ReadOutName";
43  DDSpecificsMatchesValueFilter filter{DDValue(attribute,name,0)};
44  DDFilteredView fv(cpv,filter);
45  fv.firstChild();
47  std::vector<int> temp = dbl_to_int(getDDDArray("Type",sv));
48  type_ = temp[0];
49 
50  setTimeFactor(100.);
51 
52  edm::LogInfo("FastTimerSim") << "FastTimerSD: Instantiation completed for "
53  << name << " of type " << type_;
54 }
55 
57 }
58 
59 uint32_t FastTimerSD::setDetUnitId(const G4Step * aStep) {
60 
61  //Find the depth segment
62  const G4ThreeVector& global = getGlobalEntryPoint();
63  const G4ThreeVector& local = getLocalEntryPoint();
64  int iz = (global.z() > 0) ? 1 : -1;
65  std::pair<int,int> izphi = ((ftcons) ? ((type_ == 1) ?
66  (ftcons->getZPhi(std::abs(local.z()),local.phi())) :
67  (ftcons->getEtaPhi(local.perp(),local.phi()))) :
68  (std::pair<int,int>(0,0)));
69  uint32_t id = FastTimeDetId(type_,izphi.first,izphi.second,iz).rawId();
70 #ifdef EDM_ML_DEBUG
71  edm::LogVerbatim("FastTimerSD")
72  << "Volume " << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName()
73  << ": " << global.z()
74  << " Iz(eta)phi " << izphi.first << ":" << izphi.second << ":"
75  << iz << " id " << std::hex << id << std::dec;
76 #endif
77  return id;
78 }
79 
80 void FastTimerSD::update(const BeginOfJob * job) {
81 
82  const edm::EventSetup* es = (*job)();
84  es->get<IdealGeometryRecord>().get(fdc);
85  if (fdc.isValid()) {
86  ftcons = &(*fdc);
87  } else {
88  edm::LogError("FastTimerSim") << "FastTimerSD : Cannot find FastTimeDDDConstants";
89  throw cms::Exception("Unknown", "FastTimerSD") << "Cannot find FastTimeDDDConstants\n";
90  }
91 #ifdef EDM_ML_DEBUG
92  edm::LogVerbatim("FastTimerSD")
93  << "FastTimerSD::Initialized with FastTimeDDDConstants\n";
94 #endif
95 }
96 
97 std::vector<double> FastTimerSD::getDDDArray(const std::string & str,
98  const DDsvalues_type & sv) {
99 
100  DDValue value(str);
101  if (DDfetch(&sv,value)) {
102  const std::vector<double> & fvec = value.doubles();
103  int nval = fvec.size();
104  if (nval < 1) {
105  edm::LogError("FastTimerSim") << "FastTimerSD : # of " << str
106  << " bins " << nval << " < 1 ==> illegal";
107  throw cms::Exception("DDException") << "FastTimerSD: cannot get array " << str;
108  }
109  return fvec;
110  } else {
111  edm::LogError("FastTimerSim") << "FastTimerSD: cannot get array " << str;
112  throw cms::Exception("DDException") << "FastTimerSD: cannot get array " << str;
113  }
114 }
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
FastTimerSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: FastTimerSD.cc:30
std::pair< int, int > getEtaPhi(double r, double phi) const
const G4ThreeVector & getLocalEntryPoint() const
Definition: TimingSD.h:69
#define nullptr
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
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
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: FastTimerSD.cc:80
const FastTimeDDDConstants * ftcons
Definition: FastTimerSD.h:43
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::pair< int, int > getZPhi(double z, double phi) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
uint32_t setDetUnitId(const G4Step *) override
Definition: FastTimerSD.cc:59
const G4ThreeVector & getGlobalEntryPoint() const
Definition: TimingSD.h:70
~FastTimerSD() override
Definition: FastTimerSD.cc:56
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: FastTimerSD.cc:97
DDsvalues_type mergedSpecifics() const
T get() const
Definition: EventSetup.h:71
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)