CMS 3D CMS Logo

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