CMS 3D CMS Logo

SiStripBadModuleConfigurableFakeESSource.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CalibTracker/SiStripESProducers
4 // Class: SiStripBadModuleConfigurableFakeESSource
5 //
14 // system include files
15 #include <memory>
16 
17 // user include files
20 
25 
27 public:
30 
32 
33  typedef std::unique_ptr<SiStripBadStrip> ReturnType;
34  ReturnType produce(const SiStripBadModuleRcd&);
35 
36 private:
37  using Parameters = std::vector<edm::ParameterSet>;
41  bool m_doByAPVs;
42 
43  std::vector<uint32_t> selectDetectors(const TrackerTopology* tTopo, const std::vector<uint32_t>& detIds) const;
44  std::vector<std::pair<uint32_t,std::vector<uint32_t> > > selectAPVs() const;
45 };
46 
51 
53 {
54  setWhatProduced(this);
55  findingRecord<SiStripBadModuleRcd>();
56 
57  m_badComponentList = iConfig.getUntrackedParameter<Parameters>("BadComponentList");
58  m_doByAPVs = iConfig.getUntrackedParameter<bool>("doByAPVs", false);
59  m_badAPVsList = iConfig.getUntrackedParameter<Parameters>("BadAPVList");
60  m_printDebug = iConfig.getUntrackedParameter<bool>("printDebug", false);
61 }
62 
64 
66 {
67  iValidity = edm::ValidityInterval{iov.beginOfTime(), iov.endOfTime()};
68 }
69 
70 // ------------ method called to produce the data ------------
73 {
74  using namespace edm::es;
75 
77  iRecord.getRecord<TrackerTopologyRcd>().get(tTopo);
78 
79  auto quality = std::make_unique<SiStripQuality>();
80 
82 
83  if(!m_doByAPVs){
84  std::vector<uint32_t> selDetIds{selectDetectors(tTopo.product(), reader->getAllDetIds())};
85  edm::LogInfo("SiStripQualityConfigurableFakeESSource")<<"[produce] number of selected dets to be removed " << selDetIds.size() <<std::endl;
86 
87  std::stringstream ss;
88  for ( const auto selId : selDetIds ) {
89  SiStripQuality::InputVector theSiStripVector;
90 
91  unsigned short firstBadStrip{0};
92  unsigned short NconsecutiveBadStrips = reader->getNumberOfApvsAndStripLength(selId).first * 128;
93  unsigned int theBadStripRange{quality->encode(firstBadStrip,NconsecutiveBadStrips)};
94 
95  if (m_printDebug) {
96  ss << "detid " << selId << " \t"
97  << " firstBadStrip " << firstBadStrip << "\t "
98  << " NconsecutiveBadStrips " << NconsecutiveBadStrips << "\t "
99  << " packed integer " << std::hex << theBadStripRange << std::dec
100  << std::endl;
101  }
102 
103  theSiStripVector.push_back(theBadStripRange);
104 
105  if ( ! quality->put(selId,SiStripBadStrip::Range{theSiStripVector.begin(),theSiStripVector.end()}) ) {
106  edm::LogError("SiStripQualityConfigurableFakeESSource") << "[produce] detid already exists";
107  }
108  }
109  if (m_printDebug) {
110  edm::LogInfo("SiStripQualityConfigurableFakeESSource") << ss.str();
111  }
112  quality->cleanUp();
113  //quality->fillBadComponents();
114  } else {
115 
116  std::vector<std::pair<uint32_t,std::vector<uint32_t> > > selAPVs{selectAPVs()};
117  edm::LogInfo("SiStripQualityConfigurableFakeESSource")<<"[produce] number of selected dets to be removed " << selAPVs.size() <<std::endl;
118 
119  std::stringstream ss;
120  for ( const auto selId : selAPVs ) {
121  SiStripQuality::InputVector theSiStripVector;
122  auto the_detid = selId.first;
123 
124  for( const auto apv : selId.second ) {
125 
126  unsigned short firstBadStrip = apv*128;
127  unsigned short NconsecutiveBadStrips = 128;
128  unsigned int theBadStripRange{quality->encode(firstBadStrip,NconsecutiveBadStrips)};
129 
130  if (m_printDebug) {
131  ss << "detid " << the_detid << " \t"
132  << " firstBadStrip " << firstBadStrip << "\t "
133  << " NconsecutiveBadStrips " << NconsecutiveBadStrips << "\t "
134  << " packed integer " << std::hex << theBadStripRange << std::dec
135  << std::endl;
136  }
137 
138  theSiStripVector.push_back(theBadStripRange);
139  }
140 
141  if ( ! quality->put(the_detid,SiStripBadStrip::Range{theSiStripVector.begin(),theSiStripVector.end()}) ) {
142  edm::LogError("SiStripQualityConfigurableFakeESSource") << "[produce] detid already exists";
143  }
144  } // loop on the packed list of detid/apvs
145 
146  if (m_printDebug) {
147  edm::LogInfo("SiStripQualityConfigurableFakeESSource") << ss.str();
148  }
149  quality->cleanUp();
150 
151  } // do it by APVs
152 
153  if (m_printDebug){
154  std::stringstream ss1;
155  for ( const auto& badComp : quality->getBadComponentList() ) {
156  ss1 << "bad module " << badComp.detid << " " << badComp.BadModule << "\n";
157  }
158  edm::LogInfo("SiStripQualityConfigurableFakeESSource") << ss1.str();
159  }
160 
161  return quality;
162 }
163 
164 namespace {
165  bool _isSel( uint32_t requested, uint32_t i )
166  { // internal helper: accept all i if requested is 0, otherwise require match
167  return ( requested == 0 ) || ( requested == i );
168  }
169 
170  SiStripDetId::SubDetector subDetFromString( const std::string& subDetStr )
171  {
173  if ( subDetStr=="TIB" ) subDet = SiStripDetId::TIB;
174  else if ( subDetStr=="TID" ) subDet = SiStripDetId::TID;
175  else if ( subDetStr=="TOB" ) subDet = SiStripDetId::TOB;
176  else if ( subDetStr=="TEC" ) subDet = SiStripDetId::TEC;
177  return subDet;
178  }
179 }
180 
181 std::vector<std::pair<uint32_t,std::vector<uint32_t> > > SiStripBadModuleConfigurableFakeESSource::selectAPVs() const
182 {
183 
184  std::vector<std::pair<uint32_t,std::vector<uint32_t> > > selList;
185  selList.clear();
186 
187  for ( const auto& badAPV : m_badAPVsList ) {
188  const uint32_t det{badAPV.getParameter<uint32_t>("DetId")};
189  std::vector<uint32_t> apvs{badAPV.getParameter<std::vector<uint32_t> >("APVs")};
190  auto pair = std::make_pair(det,apvs);
191  selList.push_back(pair);
192  }
193  return selList;
194 }
195 
196 
197 std::vector<uint32_t> SiStripBadModuleConfigurableFakeESSource::selectDetectors(const TrackerTopology* tTopo, const std::vector<uint32_t>& detIds) const
198 {
199  std::vector<uint32_t> selList;
200  std::stringstream ss;
201  for ( const auto& badComp : m_badComponentList ) {
202  const std::string subDetStr{badComp.getParameter<std::string>("SubDet")};
203  if (m_printDebug) ss << "Bad SubDet " << subDetStr << " \t";
204  const SiStripDetId::SubDetector subDet = subDetFromString(subDetStr);
205 
206  const std::vector<uint32_t> genericBadDetIds{badComp.getUntrackedParameter<std::vector<uint32_t>>("detidList", std::vector<uint32_t>())};
207  const bool anySubDet{ ! genericBadDetIds.empty() };
208 
209  std::cout << "genericBadDetIds.size() = " << genericBadDetIds.size() << std::endl;
210 
211  using DetIdIt = std::vector<uint32_t>::const_iterator;
212  const DetIdIt beginDetIt = std::lower_bound(detIds.begin(), detIds.end(), DetId(DetId::Tracker, anySubDet ? SiStripDetId::TIB : subDet ).rawId());
213  const DetIdIt endDetIt = std::lower_bound(detIds.begin(), detIds.end(), DetId(DetId::Tracker, anySubDet ? SiStripDetId::TEC+1 : subDet+1).rawId());
214 
215  if ( anySubDet ) {
216  std::copy_if(beginDetIt, endDetIt, std::back_inserter(selList),
217  [&genericBadDetIds] ( uint32_t detId ) {
218  std::cout << "AnySubDet" << detId << std::endl;
219  return std::find(genericBadDetIds.begin(), genericBadDetIds.end(), detId) != genericBadDetIds.end();
220  });
221  } else {
222  switch (subDet) {
223  case SiStripDetId::TIB:
224  std::copy_if(beginDetIt, endDetIt, std::back_inserter(selList),
225  [tTopo,&badComp] ( uint32_t detectorId ) {
226  const DetId detId{detectorId};
227  return ( ( detId.subdetId() == SiStripDetId::TIB )
228  && _isSel(badComp.getParameter<uint32_t>("layer") , tTopo->tibLayer(detId))
229  && _isSel(badComp.getParameter<uint32_t>("bkw_frw"), tTopo->tibIsZPlusSide(detId) ? 2 : 1 )
230  && _isSel(badComp.getParameter<uint32_t>("int_ext"), tTopo->tibIsInternalString(detId) ? 1 : 2 )
231  && _isSel(badComp.getParameter<uint32_t>("ster") , tTopo->tibIsStereo(detId) ? 1 : ( tTopo->tibIsRPhi(detId) ? 2 : -1 ))
232  && _isSel(badComp.getParameter<uint32_t>("string_"), tTopo->tibString(detId))
233  && _isSel(badComp.getParameter<uint32_t>("detid") , detId.rawId())
234  );
235  } );
236  break;
237  case SiStripDetId::TID:
238  std::copy_if(beginDetIt, endDetIt, std::back_inserter(selList),
239  [tTopo,&badComp] ( uint32_t detectorId ) {
240  const DetId detId{detectorId};
241  return ( ( detId.subdetId() == SiStripDetId::TID )
242  && _isSel(badComp.getParameter<uint32_t>("wheel"), tTopo->tidWheel(detId))
243  && _isSel(badComp.getParameter<uint32_t>("side") , tTopo->tidIsZPlusSide(detId) ? 2 : 1 )
244  && _isSel(badComp.getParameter<uint32_t>("ster") , tTopo->tidIsStereo(detId) ? 1 : ( tTopo->tidIsRPhi(detId) ? 2 : -1 ) )
245  && _isSel(badComp.getParameter<uint32_t>("ring") , tTopo->tidRing(detId))
246  && _isSel(badComp.getParameter<uint32_t>("detid"), detId.rawId())
247  );
248  } );
249  break;
250  case SiStripDetId::TOB:
251  std::copy_if(beginDetIt, endDetIt, std::back_inserter(selList),
252  [tTopo,&badComp] ( uint32_t detectorId ) {
253  const DetId detId{detectorId};
254  return ( ( detId.subdetId() == SiStripDetId::TOB )
255  && _isSel(badComp.getParameter<uint32_t>("layer") , tTopo->tobLayer(detId))
256  && _isSel(badComp.getParameter<uint32_t>("bkw_frw"), tTopo->tobIsZPlusSide(detId) ? 2 : 1)
257  && _isSel(badComp.getParameter<uint32_t>("ster") , tTopo->tobIsStereo(detId) ? 1 : ( tTopo->tobIsRPhi(detId) ? 2 : -1 ))
258  && _isSel(badComp.getParameter<uint32_t>("rod") , tTopo->tobRod(detId))
259  && _isSel(badComp.getParameter<uint32_t>("detid") , detId.rawId())
260  );
261  } );
262  break;
263  case SiStripDetId::TEC:
264  std::copy_if(beginDetIt, endDetIt, std::back_inserter(selList),
265  [tTopo,&badComp] ( uint32_t detectorId ) {
266  const DetId detId{detectorId};
267  return ( ( detId.subdetId() == SiStripDetId::TEC )
268  && _isSel(badComp.getParameter<uint32_t>("wheel") , tTopo->tecWheel(detId))
269  && _isSel(badComp.getParameter<uint32_t>("side") , tTopo->tecIsZPlusSide(detId) ? 2 : 1 )
270  && _isSel(badComp.getParameter<uint32_t>("ster") , tTopo->tecIsStereo(detId) ? 1 : 2 )
271  && _isSel(badComp.getParameter<uint32_t>("petal_bkw_frw"), tTopo->tecIsFrontPetal(detId) ? 2 : 2)
272  && _isSel(badComp.getParameter<uint32_t>("petal") , tTopo->tecPetalNumber(detId))
273  && _isSel(badComp.getParameter<uint32_t>("ring") , tTopo->tecRing(detId))
274  && _isSel(badComp.getParameter<uint32_t>("detid") , detId.rawId())
275  );
276  } );
277  break;
278  default:
279  break;
280  }
281  }
282  }
283  if (m_printDebug) {
284  edm::LogInfo("SiStripBadModuleGenerator") << ss.str();
285  }
286  return selList;
287 }
288 
289 //define this as a plug-in
auto setWhatProduced(T *iThis, const es::Label &iLabel={})
Definition: ESProducer.h:124
T getUntrackedParameter(std::string const &, T const &) const
unsigned int tibLayer(const DetId &id) const
unsigned int tibString(const DetId &id) const
unsigned int tidRing(const DetId &id) const
bool tobIsStereo(const DetId &id) const
unsigned int tecRing(const DetId &id) const
ring id
static const IOVSyncValue & endOfTime()
Definition: IOVSyncValue.cc:82
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
unsigned int tidWheel(const DetId &id) const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
bool tobIsRPhi(const DetId &id) const
bool tidIsStereo(const DetId &id) const
bool tibIsZPlusSide(const DetId &id) const
bool tecIsStereo(const DetId &id) const
bool tibIsRPhi(const DetId &id) const
static const IOVSyncValue & beginOfTime()
Definition: IOVSyncValue.cc:88
bool tobIsZPlusSide(const DetId &id) const
std::vector< std::pair< uint32_t, std::vector< uint32_t > > > selectAPVs() const
void setIntervalFor(const edm::eventsetup::EventSetupRecordKey &, const edm::IOVSyncValue &iov, edm::ValidityInterval &iValidity) override
bool tecIsFrontPetal(const DetId &id) const
Definition: DetId.h:18
bool tidIsRPhi(const DetId &id) const
bool tidIsZPlusSide(const DetId &id) const
std::vector< uint32_t > selectDetectors(const TrackerTopology *tTopo, const std::vector< uint32_t > &detIds) const
bool tibIsStereo(const DetId &id) const
bool tecIsZPlusSide(const DetId &id) const
std::pair< ContainerIterator, ContainerIterator > Range
unsigned int tecPetalNumber(const DetId &id) const
#define DEFINE_FWK_EVENTSETUP_SOURCE(type)
Definition: SourceFactory.h:91
unsigned int tobRod(const DetId &id) const
bool tibIsInternalString(const DetId &id) const
unsigned int tecWheel(const DetId &id) const
T const * product() const
Definition: ESHandle.h:86
unsigned int tobLayer(const DetId &id) const
Container InputVector