CMS 3D CMS Logo

SiStripQualityHotStripIdentifier.cc
Go to the documentation of this file.
2 
4 
12 #include <fstream>
13 #include <iostream>
14 #include <memory>
15 #include <sstream>
16 
17 //Insert here the include to the algos
19 
22  m_cacheID_(0),
23  dataLabel_(iConfig.getUntrackedParameter<std::string>("dataLabel", "")),
24  conf_(iConfig),
25  fp_(iConfig.getUntrackedParameter<edm::FileInPath>(
26  "file", edm::FileInPath("CalibTracker/SiStripCommon/data/SiStripDetInfo.dat"))),
27  Cluster_src_(iConfig.getParameter<edm::InputTag>("Cluster_src")),
28  Track_src_(iConfig.getUntrackedParameter<edm::InputTag>("Track_src")),
29  tracksCollection_in_EventTree(iConfig.getUntrackedParameter<bool>("RemoveTrackClusters", false)),
30  tTopo(nullptr) {
32 
34  MinClusterWidth_ = pset.getUntrackedParameter<uint32_t>("minWidth", 1);
35  MaxClusterWidth_ = pset.getUntrackedParameter<uint32_t>("maxWidth", 1000);
36 
37  bookHistos();
38 }
39 
41 
42 std::unique_ptr<SiStripBadStrip> SiStripQualityHotStripIdentifier::getNewObject() {
43  auto obj = std::make_unique<SiStripBadStrip>();
44 
46  std::string AlgoName = parameters.getParameter<std::string>("AlgoName");
47  if (AlgoName == "SiStripHotStripAlgorithmFromClusterOccupancy") {
48  edm::LogInfo("SiStripQualityHotStripIdentifier")
49  << " [SiStripQualityHotStripIdentifier::getNewObject] call to SiStripHotStripAlgorithmFromClusterOccupancy"
50  << std::endl;
51 
53  theIdentifier.setProbabilityThreshold(parameters.getUntrackedParameter<double>("ProbabilityThreshold", 1.E-7));
54  theIdentifier.setMinNumEntries(parameters.getUntrackedParameter<uint32_t>("MinNumEntries", 100));
55  theIdentifier.setMinNumEntriesPerStrip(parameters.getUntrackedParameter<uint32_t>("MinNumEntriesPerStrip", 5));
56 
57  SiStripQuality* qobj = new SiStripQuality();
59 
60  edm::LogInfo("SiStripQualityHotStripIdentifier")
61  << " [SiStripQualityHotStripIdentifier::getNewObject] copy SiStripObject in SiStripBadStrip" << std::endl;
62 
63  std::stringstream ss;
64 
67  for (; rIter != rIterEnd; ++rIter) {
68  SiStripBadStrip::Range range(qobj->getDataVectorBegin() + rIter->ibegin,
69  qobj->getDataVectorBegin() + rIter->iend);
70  if (!obj->put(rIter->detid, range))
71  edm::LogError("SiStripQualityHotStripIdentifier")
72  << "[SiStripQualityHotStripIdentifier::getNewObject] detid already exists" << std::endl;
73  }
74  edm::LogInfo("SiStripQualityHotStripIdentifier")
75  << " [SiStripQualityHotStripIdentifier::getNewObject] " << ss.str() << std::endl;
76 
77  } else {
78  edm::LogError("SiStripQualityHotStripIdentifier")
79  << " [SiStripQualityHotStripIdentifier::getNewObject] call for a unknow HotStrip identification algoritm"
80  << std::endl;
81 
82  std::vector<uint32_t> a;
83  SiStripBadStrip::Range range(a.begin(), a.end());
84  if (!obj->put(0xFFFFFFFF, range))
85  edm::LogError("SiStripQualityHotStripIdentifier")
86  << "[SiStripQualityHotStripIdentifier::getNewObject] detid already exists" << std::endl;
87  }
88 
89  return obj;
90 }
91 
93  //Retrieve tracker topology from geometry
95  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
96  tTopo = tTopoHandle.product();
97 
98  resetHistos();
99  unsigned long long cacheID = iSetup.get<SiStripQualityRcd>().cacheIdentifier();
100 
101  if (m_cacheID_ == cacheID)
102  return;
103 
104  m_cacheID_ = cacheID;
105 
107 }
108 
110  //Clear map
111  ClusterPositionHistoMap.clear();
112 }
113 
115  edm::LogInfo("SiStripQualityHotStripIdentifier") << " [SiStripQualityHotStripIdentifier::resetHistos] " << std::endl;
116  SiStrip::QualityHistosMap::iterator it = ClusterPositionHistoMap.begin();
117  SiStrip::QualityHistosMap::iterator iEnd = ClusterPositionHistoMap.end();
118  for (; it != iEnd; ++it) {
119  it->second->Reset();
120  }
121 }
122 
124  edm::LogInfo("SiStripQualityHotStripIdentifier") << " [SiStripQualityHotStripIdentifier::bookHistos] " << std::endl;
125  char hname[1024];
126  std::map<uint32_t, SiStripDetInfoFileReader::DetInfo>::const_iterator it = reader->getAllData().begin();
127  std::map<uint32_t, SiStripDetInfoFileReader::DetInfo>::const_iterator iEnd = reader->getAllData().end();
128  for (; it != iEnd; ++it) {
129  sprintf(hname, "h_%d", it->first);
130  SiStrip::QualityHistosMap::iterator ref = ClusterPositionHistoMap.find(it->first);
131  if (ref == ClusterPositionHistoMap.end()) {
132  ClusterPositionHistoMap[it->first] =
133  std::make_shared<TH1F>(hname, hname, it->second.nApvs * 128, -0.5, it->second.nApvs * 128 - 0.5);
134  } else
135  edm::LogError("SiStripQualityHotStripIdentifier")
136  << " [SiStripQualityHotStripIdentifier::bookHistos] DetId " << it->first
137  << " already found in map. Ignoring new data" << std::endl;
138  }
139 }
140 
142  SiStrip::QualityHistosMap::iterator ref = ClusterPositionHistoMap.find(detid);
143  if (ref != ClusterPositionHistoMap.end())
144  ref->second->Fill(value);
145  else
146  edm::LogError("SiStripQualityHotStripIdentifier")
147  << " [SiStripQualityHotStripIdentifier::fillHisto] Histogram not found in the list for DetId " << detid
148  << " Ignoring data value " << value << std::endl;
149 }
150 
153  e.getByLabel(Cluster_src_, dsv_SiStripCluster);
154 
157  e.getByLabel(Track_src_, trackCollection);
158  if (!trackCollection.isValid()) {
159  edm::LogError("SiStripQualityHotStripIdentifier")
160  << " [SiStripQualityHotStripIdentifier::algoAnalyze] missing trackCollection with label " << Track_src_
161  << std::endl;
162  }
163  }
164 
165  std::set<const void*> vPSiStripCluster;
166  //Perform track study
168  const reco::TrackCollection tC = *(trackCollection.product());
169  int i = 0;
170  for (reco::TrackCollection::const_iterator track = tC.begin(); track != tC.end(); track++) {
171  LogTrace("SiStripQualityHotStripIdentifier")
172  << "Track number " << i + 1 << "\n\tmomentum: " << track->momentum() << "\n\tPT: " << track->pt()
173  << "\n\tvertex: " << track->vertex() << "\n\timpact parameter: " << track->d0()
174  << "\n\tcharge: " << track->charge() << "\n\tnormalizedChi2: " << track->normalizedChi2()
175  << "\n\tFrom EXTRA : "
176  << "\n\t\touter PT " << track->outerPt() << std::endl;
177 
178  //Loop on rechits
179  for (auto const& recHit : track->recHits()) {
180  if (!recHit->isValid()) {
181  LogTrace("SiStripQualityHotStripIdentifier") << "\t\t Invalid Hit " << std::endl;
182  continue;
183  }
184 
185  const SiStripRecHit2D* singleHit = dynamic_cast<const SiStripRecHit2D*>(recHit);
186  const SiStripMatchedRecHit2D* matchedHit = dynamic_cast<const SiStripMatchedRecHit2D*>(recHit);
187  const ProjectedSiStripRecHit2D* projectedHit = dynamic_cast<const ProjectedSiStripRecHit2D*>(recHit);
188 
189  if (matchedHit) {
190  vPSiStripCluster.insert((void*)&(matchedHit->monoCluster()));
191  vPSiStripCluster.insert((void*)&(matchedHit->stereoCluster()));
192  } else if (projectedHit) {
193  vPSiStripCluster.insert((void*)&*(projectedHit->originalHit().cluster()));
194  } else if (singleHit) {
195  vPSiStripCluster.insert((void*)&*(singleHit->cluster()));
196  } else {
197  LogTrace("SiStripQualityHotStripIdentifier") << "NULL hit" << std::endl;
198  }
199  }
200  }
201  }
202 
203  std::stringstream ss;
204  //Loop on Det Clusters
205  edm::DetSetVector<SiStripCluster>::const_iterator DSViter = dsv_SiStripCluster->begin();
206  for (; DSViter != dsv_SiStripCluster->end(); DSViter++) {
208  edm::DetSet<SiStripCluster>::const_iterator ClusIterEnd = DSViter->data.end();
209  for (; ClusIter != ClusIterEnd; ++ClusIter) {
210  if (MinClusterWidth_ <= ClusIter->amplitudes().size() && ClusIter->amplitudes().size() <= MaxClusterWidth_) {
211  if (std::find(vPSiStripCluster.begin(), vPSiStripCluster.end(), (void*)&*ClusIter) == vPSiStripCluster.end()) {
212  if (edm::isDebugEnabled())
213  ss << " adding cluster to histo for detid " << DSViter->id << " with barycenter " << ClusIter->barycenter()
214  << std::endl;
215  fillHisto(DSViter->id, ClusIter->barycenter());
216  }
217  }
218  }
219  }
220  LogTrace("SiStripQualityHotStripIdentifier") << ss.str();
221 }
size
Write out results.
T getParameter(std::string const &) const
bool isDebugEnabled()
T getUntrackedParameter(std::string const &, T const &) const
void fillHisto(uint32_t detid, float value)
void extractBadStrips(SiStripQuality *, HistoMap &, edm::ESHandle< SiStripQuality > &)
SiStripCluster const & monoCluster() const
#define nullptr
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
Registry::const_iterator RegistryIterator
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::unique_ptr< SiStripBadStrip > getNewObject() override
void algoBeginRun(const edm::Run &, const edm::EventSetup &) override
size_type size() const
Definition: DetSet.h:62
SiStripQualityHotStripIdentifier(const edm::ParameterSet &)
edmNew::DetSet< SiStripCluster >::const_iterator ClusIter
RegistryIterator getRegistryVectorEnd() const
void algoAnalyze(const edm::Event &, const edm::EventSetup &) override
ClusterRef cluster() const
edm::ESHandle< SiStripQuality > SiStripQuality_
Definition: value.py:1
bool isValid() const
Definition: HandleBase.h:70
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:488
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:325
#define LogTrace(id)
ContainerIterator getDataVectorBegin() const
SiStripRecHit2D originalHit() const
T const * product() const
Definition: Handle.h:69
RegistryIterator getRegistryVectorBegin() const
SiStripCluster const & stereoCluster() const
HLT enums.
double a
Definition: hdecay.h:119
std::pair< ContainerIterator, ContainerIterator > Range
T get() const
Definition: EventSetup.h:73
std::string fullPath() const
Definition: FileInPath.cc:163
const std::map< uint32_t, DetInfo > & getAllData() const
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:314
collection_type::const_iterator const_iterator
Definition: DetSet.h:32
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:102
T const * product() const
Definition: ESHandle.h:86
Definition: Run.h:45