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::shared_ptr<SiStripBadStrip> ReturnType;
34  ReturnType produce(const SiStripBadModuleRcd&);
35 
36 private:
37  using Parameters = std::vector<edm::ParameterSet>;
41 
42  std::vector<uint32_t> selectDetectors(const TrackerTopology* tTopo, const std::vector<uint32_t>& detIds) const;
43 };
44 
49 
51 {
52  setWhatProduced(this);
53  findingRecord<SiStripBadModuleRcd>();
54 
55  m_badComponentList = iConfig.getUntrackedParameter<Parameters>("BadComponentList");
56  m_file = iConfig.getParameter<edm::FileInPath>("file");
57  m_printDebug = iConfig.getUntrackedParameter<bool>("printDebug", false);
58 }
59 
61 
63 {
64  iValidity = edm::ValidityInterval{iov.beginOfTime(), iov.endOfTime()};
65 }
66 
67 // ------------ method called to produce the data ------------
70 {
71  using namespace edm::es;
72 
74  iRecord.getRecord<TrackerTopologyRcd>().get(tTopo);
75 
76  std::shared_ptr<SiStripQuality> quality{new SiStripQuality};
77 
79 
80  std::vector<uint32_t> selDetIds{selectDetectors(tTopo.product(), reader.getAllDetIds())};
81  edm::LogInfo("SiStripQualityConfigurableFakeESSource")<<"[produce] number of selected dets to be removed " << selDetIds.size() <<std::endl;
82 
83  std::stringstream ss;
84  for ( const auto selId : selDetIds ) {
85  SiStripQuality::InputVector theSiStripVector;
86 
87  unsigned short firstBadStrip{0};
88  unsigned short NconsecutiveBadStrips = reader.getNumberOfApvsAndStripLength(selId).first * 128;
89  unsigned int theBadStripRange{quality->encode(firstBadStrip,NconsecutiveBadStrips)};
90 
91  if (m_printDebug) {
92  ss << "detid " << selId << " \t"
93  << " firstBadStrip " << firstBadStrip << "\t "
94  << " NconsecutiveBadStrips " << NconsecutiveBadStrips << "\t "
95  << " packed integer " << std::hex << theBadStripRange << std::dec
96  << std::endl;
97  }
98 
99  theSiStripVector.push_back(theBadStripRange);
100 
101  if ( ! quality->put(selId,SiStripBadStrip::Range{theSiStripVector.begin(),theSiStripVector.end()}) ) {
102  edm::LogError("SiStripQualityConfigurableFakeESSource") << "[produce] detid already exists";
103  }
104  }
105  if (m_printDebug) {
106  edm::LogInfo("SiStripQualityConfigurableFakeESSource") << ss.str();
107  }
108  quality->cleanUp();
109  //quality->fillBadComponents();
110 
111  if (m_printDebug){
112  std::stringstream ss1;
113  for ( const auto& badComp : quality->getBadComponentList() ) {
114  ss1 << "bad module " << badComp.detid << " " << badComp.BadModule << "\n";
115  }
116  edm::LogInfo("SiStripQualityConfigurableFakeESSource") << ss1.str();
117  }
118 
119  return quality;
120 }
121 
122 namespace {
123  bool _isSel( uint32_t requested, uint32_t i )
124  { // internal helper: accept all i if requested is 0, otherwise require match
125  return ( requested == 0 ) || ( requested == i );
126  }
127 
128  SiStripDetId::SubDetector subDetFromString( const std::string& subDetStr )
129  {
131  if ( subDetStr=="TIB" ) subDet = SiStripDetId::TIB;
132  else if ( subDetStr=="TID" ) subDet = SiStripDetId::TID;
133  else if ( subDetStr=="TOB" ) subDet = SiStripDetId::TOB;
134  else if ( subDetStr=="TEC" ) subDet = SiStripDetId::TEC;
135  return subDet;
136  }
137 }
138 
139 std::vector<uint32_t> SiStripBadModuleConfigurableFakeESSource::selectDetectors(const TrackerTopology* tTopo, const std::vector<uint32_t>& detIds) const
140 {
141  std::vector<uint32_t> selList;
142  std::stringstream ss;
143  for ( const auto& badComp : m_badComponentList ) {
144  const std::string subDetStr{badComp.getParameter<std::string>("SubDet")};
145  if (m_printDebug) ss << "Bad SubDet " << subDetStr << " \t";
146  const SiStripDetId::SubDetector subDet = subDetFromString(subDetStr);
147 
148  const std::vector<uint32_t> genericBadDetIds{badComp.getUntrackedParameter<std::vector<uint32_t>>("detidList", std::vector<uint32_t>())};
149  const bool anySubDet{ ! genericBadDetIds.empty() };
150 
151  std::cout << "genericBadDetIds.size() = " << genericBadDetIds.size() << std::endl;
152 
153  using DetIdIt = std::vector<uint32_t>::const_iterator;
154  const DetIdIt beginDetIt = std::lower_bound(detIds.begin(), detIds.end(), DetId(DetId::Tracker, anySubDet ? SiStripDetId::TIB : subDet ).rawId());
155  const DetIdIt endDetIt = std::lower_bound(detIds.begin(), detIds.end(), DetId(DetId::Tracker, anySubDet ? SiStripDetId::TEC+1 : subDet+1).rawId());
156 
157  if ( anySubDet ) {
158  std::copy_if(beginDetIt, endDetIt, std::back_inserter(selList),
159  [&genericBadDetIds] ( uint32_t detId ) {
160  std::cout << "AnySubDet" << detId << std::endl;
161  return std::find(genericBadDetIds.begin(), genericBadDetIds.end(), detId) != genericBadDetIds.end();
162  });
163  } else {
164  switch (subDet) {
165  case SiStripDetId::TIB:
166  std::copy_if(beginDetIt, endDetIt, std::back_inserter(selList),
167  [tTopo,&badComp] ( uint32_t detectorId ) {
168  const DetId detId{detectorId};
169  return ( ( detId.subdetId() == SiStripDetId::TIB )
170  && _isSel(badComp.getParameter<uint32_t>("layer") , tTopo->tibLayer(detId))
171  && _isSel(badComp.getParameter<uint32_t>("bkw_frw"), tTopo->tibIsZPlusSide(detId) ? 2 : 1 )
172  && _isSel(badComp.getParameter<uint32_t>("int_ext"), tTopo->tibIsInternalString(detId) ? 1 : 2 )
173  && _isSel(badComp.getParameter<uint32_t>("ster") , tTopo->tibIsStereo(detId) ? 1 : ( tTopo->tibIsRPhi(detId) ? 2 : -1 ))
174  && _isSel(badComp.getParameter<uint32_t>("string_"), tTopo->tibString(detId))
175  && _isSel(badComp.getParameter<uint32_t>("detid") , detId.rawId())
176  );
177  } );
178  break;
179  case SiStripDetId::TID:
180  std::copy_if(beginDetIt, endDetIt, std::back_inserter(selList),
181  [tTopo,&badComp] ( uint32_t detectorId ) {
182  const DetId detId{detectorId};
183  return ( ( detId.subdetId() == SiStripDetId::TID )
184  && _isSel(badComp.getParameter<uint32_t>("wheel"), tTopo->tidWheel(detId))
185  && _isSel(badComp.getParameter<uint32_t>("side") , tTopo->tidIsZPlusSide(detId) ? 2 : 1 )
186  && _isSel(badComp.getParameter<uint32_t>("ster") , tTopo->tidIsStereo(detId) ? 1 : ( tTopo->tidIsRPhi(detId) ? 2 : -1 ) )
187  && _isSel(badComp.getParameter<uint32_t>("ring") , tTopo->tidRing(detId))
188  && _isSel(badComp.getParameter<uint32_t>("detid"), detId.rawId())
189  );
190  } );
191  break;
192  case SiStripDetId::TOB:
193  std::copy_if(beginDetIt, endDetIt, std::back_inserter(selList),
194  [tTopo,&badComp] ( uint32_t detectorId ) {
195  const DetId detId{detectorId};
196  return ( ( detId.subdetId() == SiStripDetId::TOB )
197  && _isSel(badComp.getParameter<uint32_t>("layer") , tTopo->tobLayer(detId))
198  && _isSel(badComp.getParameter<uint32_t>("bkw_frw"), tTopo->tobIsZPlusSide(detId) ? 2 : 1)
199  && _isSel(badComp.getParameter<uint32_t>("ster") , tTopo->tobIsStereo(detId) ? 1 : ( tTopo->tobIsRPhi(detId) ? 2 : -1 ))
200  && _isSel(badComp.getParameter<uint32_t>("rod") , tTopo->tobRod(detId))
201  && _isSel(badComp.getParameter<uint32_t>("detid") , detId.rawId())
202  );
203  } );
204  break;
205  case SiStripDetId::TEC:
206  std::copy_if(beginDetIt, endDetIt, std::back_inserter(selList),
207  [tTopo,&badComp] ( uint32_t detectorId ) {
208  const DetId detId{detectorId};
209  return ( ( detId.subdetId() == SiStripDetId::TEC )
210  && _isSel(badComp.getParameter<uint32_t>("wheel") , tTopo->tecWheel(detId))
211  && _isSel(badComp.getParameter<uint32_t>("side") , tTopo->tecIsZPlusSide(detId) ? 2 : 1 )
212  && _isSel(badComp.getParameter<uint32_t>("ster") , tTopo->tecIsStereo(detId) ? 1 : 2 )
213  && _isSel(badComp.getParameter<uint32_t>("petal_bkw_frw"), tTopo->tecIsFrontPetal(detId) ? 2 : 2)
214  && _isSel(badComp.getParameter<uint32_t>("petal") , tTopo->tecPetalNumber(detId))
215  && _isSel(badComp.getParameter<uint32_t>("ring") , tTopo->tecRing(detId))
216  && _isSel(badComp.getParameter<uint32_t>("detid") , detId.rawId())
217  );
218  } );
219  break;
220  default:
221  break;
222  }
223  }
224  }
225  if (m_printDebug) {
226  edm::LogInfo("SiStripBadModuleGenerator") << ss.str();
227  }
228  return selList;
229 }
230 
231 //define this as a plug-in
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
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:97
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:7
bool tobIsRPhi(const DetId &id) const
void setWhatProduced(T *iThis, const es::Label &iLabel=es::Label())
Definition: ESProducer.h:115
bool tidIsStereo(const DetId &id) const
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
bool tibIsZPlusSide(const DetId &id) const
bool tecIsStereo(const DetId &id) const
void setIntervalFor(const edm::eventsetup::EventSetupRecordKey &, const edm::IOVSyncValue &iov, edm::ValidityInterval &iValidity)
bool tibIsRPhi(const DetId &id) const
static const IOVSyncValue & beginOfTime()
bool tobIsZPlusSide(const DetId &id) const
bool tecIsFrontPetal(const DetId &id) const
#define DEFINE_FWK_EVENTSETUP_SOURCE(type)
Definition: SourceFactory.h:92
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
std::string fullPath() const
Definition: FileInPath.cc:184
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