CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiPixelBadModuleByHandBuilder.cc
Go to the documentation of this file.
1 // system includes
2 #include <vector>
3 #include <cmath>
4 #include <fstream>
5 #include <iostream>
6 
7 // user includes
21 
22 class SiPixelBadModuleByHandBuilder : public ConditionDBWriter<SiPixelQuality> {
23 public:
26 
27 private:
28  std::unique_ptr<SiPixelQuality> getNewObject() override;
29 
30  void algoBeginRun(const edm::Run& run, const edm::EventSetup& es) override {
31  if (!tTopo_) {
32  tTopo_ = std::make_unique<TrackerTopology>(es.getData(tkTopoToken_));
33  }
34  };
35 
36 private:
38  const bool printdebug_;
39  typedef std::vector<edm::ParameterSet> Parameters;
42  std::unique_ptr<TrackerTopology> tTopo_;
43 };
44 
47  tkTopoToken_(esConsumes<edm::Transition::BeginRun>()),
48  printdebug_(iConfig.getUntrackedParameter<bool>("printDebug", false)),
49  BadModuleList_(iConfig.getUntrackedParameter<Parameters>("BadModuleList")),
50  ROCListFile_(iConfig.getUntrackedParameter<std::string>("ROCListFile")) {}
51 
53 
54 std::unique_ptr<SiPixelQuality> SiPixelBadModuleByHandBuilder::getNewObject() {
55  auto obj = std::make_unique<SiPixelQuality>();
56 
57  for (Parameters::iterator it = BadModuleList_.begin(); it != BadModuleList_.end(); ++it) {
58  if (printdebug_) {
59  edm::LogInfo("SiPixelBadModuleByHandBuilder") << " BadModule " << *it << " \t" << std::endl;
60  }
61 
63  BadModule.errorType = 3;
64  BadModule.BadRocs = 0;
65  BadModule.DetID = it->getParameter<uint32_t>("detid");
66  std::string errorstring = it->getParameter<std::string>("errortype");
67  if (printdebug_) {
68  edm::LogInfo("SiPixelBadModuleByHandBuilder")
69  << "now looking at detid " << BadModule.DetID << ", string " << errorstring << std::endl;
70  }
71 
73  // errortype "whole" = int 0 in DB //
74  // errortype "tbmA" = int 1 in DB //
75  // errortype "tbmB" = int 2 in DB //
76  // errortype "none" = int 3 in DB //
78 
80  //each bad roc correspond to a bit to 1: num= //
81  // 0 <-> all good rocs //
82  // 1 <-> only roc 0 bad //
83  // 2<-> only roc 1 bad //
84  // 3<-> roc 0 and 1 bad //
85  // 4 <-> only roc 2 bad //
86  // ... //
88 
89  if (errorstring == "whole") {
90  BadModule.errorType = 0;
91  BadModule.BadRocs = 65535;
92  } //corresponds to all rocs being bad
93  else if (errorstring == "tbmA") {
94  BadModule.errorType = 1;
95  BadModule.BadRocs = 255;
96  } //corresponds to Rocs 0-7 being bad
97  else if (errorstring == "tbmB") {
98  BadModule.errorType = 2;
99  BadModule.BadRocs = 65280;
100  } //corresponds to Rocs 8-15 being bad
101  else if (errorstring == "none") {
102  BadModule.errorType = 3;
103  // badroclist_ = iConfig.getUntrackedParameter<std::vector<uint32_t> >("badroclist");
104  std::vector<uint32_t> BadRocList = it->getParameter<std::vector<uint32_t> >("badroclist");
105  short badrocs = 0;
106  for (std::vector<uint32_t>::iterator iter = BadRocList.begin(); iter != BadRocList.end(); ++iter) {
107  badrocs += 1 << *iter; // 1 << *iter = 2^{*iter} using bitwise shift
108  }
109  BadModule.BadRocs = badrocs;
110  }
111 
112  else
113  edm::LogError("SiPixelQuality") << "trying to fill error type " << errorstring << ", which is not defined!";
114  obj->addDisabledModule(BadModule);
115  }
116 
117  // fill DB from DQM list
118  if (!ROCListFile_.empty()) {
119  std::map<uint32_t, uint32_t> disabledModules;
120  std::ifstream aFile(ROCListFile_.c_str());
121  std::string aLine;
122  while (std::getline(aFile, aLine)) {
123  char name[100];
124  int roc;
125  sscanf(aLine.c_str(), "%s %d", name, &roc);
126  uint32_t detId;
127  if (name[0] == 'B') {
128  PixelBarrelName bn(name, true);
129  detId = bn.getDetId(tTopo_.get());
130  } else {
131  PixelEndcapName en(name, true);
132  detId = en.getDetId(tTopo_.get());
133  }
134  std::map<uint32_t, uint32_t>::iterator it = disabledModules.find(detId);
135  if (it == disabledModules.end())
136  it = disabledModules.insert(disabledModules.begin(), std::make_pair(detId, 0));
137  it->second |= 1 << roc;
138  //edm::LogPrint("SiPixelBadModuleByHandBuilder")<<"New module read "<<name<<" "<<roc<<" --> "<<detId<<" "<<std::bitset<32>(it->second)<<std::endl;
139  }
140 
141  for (const auto& it : disabledModules) {
143  BadModule.DetID = it.first;
144  if (it.second == 65535) { // "whole"
145  BadModule.errorType = 0;
146  BadModule.BadRocs = 65535;
147  } //corresponds to all rocs being bad
148  else if (it.second == 255) { // "tbmA"
149  BadModule.errorType = 1;
150  BadModule.BadRocs = 255;
151  } //corresponds to Rocs 0-7 being bad
152  else if (it.second == 65280) { // "tbmB"
153  BadModule.errorType = 2;
154  BadModule.BadRocs = 65280;
155  } //corresponds to Rocs 8-15 being bad
156  else { // "none"
157  BadModule.errorType = 3;
158  BadModule.BadRocs = it.second;
159  }
160 
161  obj->addDisabledModule(BadModule);
162  if (printdebug_) {
163  edm::LogVerbatim("SiPixelBadModuleByHandBuilder")
164  << "New module added: " << tTopo_->print(BadModule.DetID) << ", errorType: " << BadModule.errorType
165  << ", BadRocs: " << std::bitset<16>(it.second) << std::endl;
166  }
167  }
168  }
169 
170  return obj;
171 }
Log< level::Info, true > LogVerbatim
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void algoBeginRun(const edm::Run &run, const edm::EventSetup &es) override
SiPixelBadModuleByHandBuilder(const edm::ParameterSet &)
Log< level::Error, false > LogError
std::vector< edm::ParameterSet > Parameters
bool getData(T &iHolder) const
Definition: EventSetup.h:128
Transition
Definition: Transition.h:12
~SiPixelBadModuleByHandBuilder() override
PXFDetId getDetId()
return DetId
Log< level::Info, false > LogInfo
std::unique_ptr< SiPixelQuality > getNewObject() override
std::unique_ptr< TrackerTopology > tTopo_
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tkTopoToken_
PXBDetId getDetId()
return the DetId
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
Definition: Run.h:45