CMS 3D CMS Logo

SiStripBadComponentsDQMServiceReader.cc
Go to the documentation of this file.
3 
8 
9 #include <iostream>
10 #include <cstdio>
11 #include <sys/time.h>
12 
13 #include <boost/lexical_cast.hpp>
14 
15 using namespace std;
16 
18  : printdebug_(iConfig.getUntrackedParameter<bool>("printDebug", true)) {}
19 
21 
23  //Retrieve tracker topology from geometry
25  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
26  const TrackerTopology* const tTopo = tTopoHandle.product();
27 
28  uint32_t FedErrorMask = 1; // bit 0
29  uint32_t DigiErrorMask = 2; // bit 1
30  uint32_t ClusterErrorMask = 4; // bit 2
31 
32  edm::ESHandle<SiStripBadStrip> SiStripBadStrip_;
33  iSetup.get<SiStripBadStripRcd>().get(SiStripBadStrip_);
34  edm::LogInfo("SiStripBadComponentsDQMServiceReader")
35  << "[SiStripBadComponentsDQMServiceReader::analyze] End Reading SiStripBadStrip" << std::endl;
36 
37  std::vector<uint32_t> detid;
38  SiStripBadStrip_->getDetIds(detid);
39 
40  std::stringstream ss;
41 
42  // ss << " detid" << " \t\t\t" << "FED error" << " \t" << "Digi test failed" << " \t" << "Cluster test failed" << std::endl;
43 
44  ss << "subdet layer stereo side \t detId \t\t Errors" << std::endl;
45 
46  for (size_t id = 0; id < detid.size(); id++) {
47  SiStripBadStrip::Range range = SiStripBadStrip_->getRange(detid[id]);
48 
49  for (int it = 0; it < range.second - range.first; it++) {
50  unsigned int value = (*(range.first + it));
51  ss << detIdToString(detid[id], tTopo) << "\t" << detid[id] << "\t";
52 
53  uint32_t flag = boost::lexical_cast<uint32_t>(SiStripBadStrip_->decode(value).flag);
54 
55  printError(ss, ((flag & FedErrorMask) == FedErrorMask), "Fed error, ");
56  printError(ss, ((flag & DigiErrorMask) == DigiErrorMask), "Digi error, ");
57  printError(ss, ((flag & ClusterErrorMask) == ClusterErrorMask), "Cluster error");
58  ss << std::endl;
59 
60  if (printdebug_) {
61  ss << " firstBadStrip " << SiStripBadStrip_->decode(value).firstStrip << "\t "
62  << " NconsecutiveBadStrips " << SiStripBadStrip_->decode(value).range << "\t " // << std::endl;
63  << " flag " << SiStripBadStrip_->decode(value).flag << "\t "
64  << " packed integer " << std::hex << value << std::dec << "\t " << std::endl;
65  }
66  }
67  ss << std::endl;
68  }
69  edm::LogInfo("SiStripBadComponentsDQMServiceReader") << ss.str();
70 }
71 
73  const bool error,
74  const std::string& errorText) {
75  if (error) {
76  ss << errorText << "\t ";
77  } else {
78  ss << "\t\t ";
79  }
80 }
81 
84  int layer = 0;
85  int stereo = 0;
86  int side = -1;
87 
88  // Using the operator[] if the element does not exist it is created with the default value. That is 0 for integral types.
89  switch (detid.subdetId()) {
90  case StripSubdetector::TIB: {
91  detector = "TIB";
92  layer = tTopo->tibLayer(detid.rawId());
93  stereo = tTopo->tibStereo(detid.rawId());
94  break;
95  }
96  case StripSubdetector::TOB: {
97  detector = "TOB";
98  layer = tTopo->tobLayer(detid.rawId());
99  stereo = tTopo->tobStereo(detid.rawId());
100  break;
101  }
102  case StripSubdetector::TEC: {
103  // is this module in TEC+ or TEC-?
104  side = tTopo->tecSide(detid.rawId());
105  detector = "TEC";
106  layer = tTopo->tecWheel(detid.rawId());
107  stereo = tTopo->tecStereo(detid.rawId());
108  break;
109  }
110  case StripSubdetector::TID: {
111  // is this module in TID+ or TID-?
112  side = tTopo->tidSide(detid.rawId());
113  detector = "TID";
114  layer = tTopo->tidWheel(detid.rawId());
115  stereo = tTopo->tidStereo(detid.rawId());
116  break;
117  }
118  }
119  std::string name(detector + "\t" + std::to_string(layer) + "\t" + std::to_string(stereo) + "\t");
120  if (side == 1) {
121  name += "-";
122  } else if (side == 2) {
123  name += "+";
124  }
125  // if( side != -1 ) {
126  // name += boost::lexical_cast<string>(side);
127  // }
128 
129  return name;
130 }
unsigned short range
static constexpr auto TEC
void getDetIds(std::vector< uint32_t > &DetIds_) const
unsigned int tibLayer(const DetId &id) const
void printError(std::stringstream &ss, const bool error, const std::string &errorText)
uint32_t tobStereo(const DetId &id) const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
unsigned int tidWheel(const DetId &id) const
unsigned int tidSide(const DetId &id) const
uint32_t tidStereo(const DetId &id) const
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
Definition: value.py:1
static constexpr auto TOB
Definition: DetId.h:17
static constexpr auto TIB
unsigned short firstStrip
std::string detIdToString(const DetId &detid, const TrackerTopology *tTopo)
const Range getRange(const uint32_t detID) const
std::pair< ContainerIterator, ContainerIterator > Range
T get() const
Definition: EventSetup.h:73
uint32_t tecStereo(const DetId &id) const
void analyze(const edm::Event &, const edm::EventSetup &) override
uint32_t tibStereo(const DetId &id) const
unsigned int tecWheel(const DetId &id) const
static constexpr auto TID
T const * product() const
Definition: ESHandle.h:86
data decode(const unsigned int &value) const
unsigned int tobLayer(const DetId &id) const
unsigned int tecSide(const DetId &id) const