CMS 3D CMS Logo

SiStripBaselineAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiStripBaselineAnalyzer
4 // Class: SiStripBaselineAnalyzer
5 //
13 //
14 // Original Author: Ivan Amos Cali
15 // Created: Mon Jul 28 14:10:52 CEST 2008
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 #include <iostream>
23 
24 // user include files
28 
36 
41 
44 
48 
52 
54 
55 //ROOT inclusion
56 #include "TROOT.h"
57 #include "TFile.h"
58 #include "TNtuple.h"
59 #include "TMath.h"
60 #include "TCanvas.h"
61 #include "TH1F.h"
62 #include "TH2F.h"
63 #include "TProfile.h"
64 #include "TList.h"
65 #include "TString.h"
66 #include "TStyle.h"
67 #include "TGraph.h"
68 #include "TMultiGraph.h"
69 #include "THStack.h"
70 
71 
72 //
73 // class decleration
74 //
75 
76 class SiStripBaselineAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
77  public:
80 
81 
82  private:
83  virtual void beginJob() override ;
84  virtual void analyze(const edm::Event&, const edm::EventSetup&) override;
85  virtual void endJob() override ;
86 
87  std::auto_ptr<SiStripPedestalsSubtractor> subtractorPed_;
89  std::vector<int> pedestals;
90  uint32_t peds_cache_id;
91 
96  bool plotAPVCM_;
98 
103 
105 
107 
109  TH1F* h1Baseline_;
110  TH1F* h1Clusters_;
111  TH1F* h1APVCM_;
112  TH1F* h1Pedestals_;
113 
114  TCanvas* Canvas_;
115  std::vector<TH1F> vProcessedRawDigiHisto_;
116  std::vector<TH1F> vBaselineHisto_;
117  std::vector<TH1F> vBaselinePointsHisto_;
118  std::vector<TH1F> vClusterHisto_;
119 
121  uint16_t actualModule_;
122 };
123 
124 
126  usesResource("TFileService");
127 
128  srcBaseline_ = conf.getParameter<edm::InputTag>( "srcBaseline" );
129  srcBaselinePoints_ = conf.getParameter<edm::InputTag>( "srcBaselinePoints" );
130  srcProcessedRawDigi_ = conf.getParameter<edm::InputTag>( "srcProcessedRawDigi" );
131  srcAPVCM_ = conf.getParameter<edm::InputTag>( "srcAPVCM" );
133  nModuletoDisplay_ = conf.getParameter<uint32_t>( "nModuletoDisplay" );
134  plotClusters_ = conf.getParameter<bool>( "plotClusters" );
135  plotBaseline_ = conf.getParameter<bool>( "plotBaseline" );
136  plotBaselinePoints_ = conf.getParameter<bool>( "plotBaselinePoints" );
137  plotRawDigi_ = conf.getParameter<bool>( "plotRawDigi" );
138  plotAPVCM_ = conf.getParameter<bool>( "plotAPVCM" );
139  plotPedestals_ = conf.getParameter<bool>( "plotPedestals" );
140 
141  h1BadAPVperEvent_ = fs_->make<TH1F>("BadAPV/Event","BadAPV/Event", 2001, -0.5, 2000.5);
142  h1BadAPVperEvent_->SetXTitle("# Modules with Bad APVs");
143  h1BadAPVperEvent_->SetYTitle("Entries");
144  h1BadAPVperEvent_->SetLineWidth(2);
145  h1BadAPVperEvent_->SetLineStyle(2);
146 
147  h1APVCM_ = fs_->make<TH1F>("APV CM","APV CM", 2048, -1023.5, 1023.5);
148  h1APVCM_->SetXTitle("APV CM [adc]");
149  h1APVCM_->SetYTitle("Entries");
150  h1APVCM_->SetLineWidth(2);
151  h1APVCM_->SetLineStyle(2);
152 
153  h1Pedestals_ = fs_->make<TH1F>("Pedestals","Pedestals", 2048, -1023.5, 1023.5);
154  h1Pedestals_->SetXTitle("Pedestals [adc]");
155  h1Pedestals_->SetYTitle("Entries");
156  h1Pedestals_->SetLineWidth(2);
157  h1Pedestals_->SetLineStyle(2);
158 
159 
160 
161 }
162 
163 
165 {
166 
167 
168 
169 }
170 
171 void
173 {
174  using namespace edm;
175  if(plotPedestals_&&actualModule_ ==0){
176  uint32_t p_cache_id = es.get<SiStripPedestalsRcd>().cacheIdentifier();
177  if(p_cache_id != peds_cache_id) {
179  peds_cache_id = p_cache_id;
180  }
181 
182 
183  std::vector<uint32_t> detIdV;
184  pedestalsHandle->getDetIds(detIdV);
185 
186  for(uint32_t i=0; i < detIdV.size(); ++i){
187  pedestals.clear();
188  SiStripPedestals::Range pedestalsRange = pedestalsHandle->getRange(detIdV[i]);
189  pedestals.resize((pedestalsRange.second- pedestalsRange.first)*8/10);
190  pedestalsHandle->allPeds(pedestals, pedestalsRange);
191  for(uint32_t it=0; it < pedestals.size(); ++it) h1Pedestals_->Fill(pedestals[it]);
192  }
193  }
194 
195  if(plotAPVCM_){
197  edm::InputTag CMLabel("siStripZeroSuppression:APVCM");
198  e.getByLabel(srcAPVCM_,moduleCM);
199 
200  edm::DetSetVector<SiStripProcessedRawDigi>::const_iterator itCMDetSetV =moduleCM->begin();
201  for (; itCMDetSetV != moduleCM->end(); ++itCMDetSetV){
203  for(;itCM != itCMDetSetV->end(); ++itCM) h1APVCM_->Fill(itCM->adc());
204  }
205  }
206 
207  if(!plotRawDigi_) return;
208  subtractorPed_->init(es);
209 
210 
211 
213  e.getByLabel(srcProcessedRawDigi_,moduleRawDigi);
214 
216  if(plotBaseline_) e.getByLabel(srcBaseline_, moduleBaseline);
217 
218  edm::Handle<edm::DetSetVector<SiStripDigi> > moduleBaselinePoints;
219  if(plotBaselinePoints_) e.getByLabel(srcBaseline_, moduleBaselinePoints);
220 
222  if(plotClusters_){
223  edm::InputTag clusLabel("siStripClusters");
224  e.getByLabel(clusLabel, clusters);
225  }
226 
227  char detIds[20];
228  char evs[20];
229  char runs[20];
230 
231 
232  TFileDirectory sdProcessedRawDigis_= fs_->mkdir("ProcessedRawDigis");
233  TFileDirectory sdBaseline_= fs_->mkdir("Baseline");
234  TFileDirectory sdBaselinePoints_= fs_->mkdir("BaselinePoints");
235  TFileDirectory sdClusters_= fs_->mkdir("Clusters");
236 
237 
239  if(plotBaseline_) itDSBaseline = moduleBaseline->begin();
240  edm::DetSetVector<SiStripRawDigi>::const_iterator itRawDigis = moduleRawDigi->begin();
241 
242  uint32_t NBabAPVs = moduleRawDigi->size();
243  std::cout<< "Number of module with HIP in this event: " << NBabAPVs << std::endl;
244  h1BadAPVperEvent_->Fill(NBabAPVs);
245 
246  for (; itRawDigis != moduleRawDigi->end(); ++itRawDigis) {
247  if(actualModule_ > nModuletoDisplay_) return;
248  uint32_t detId = itRawDigis->id;
249 
250  if(plotBaseline_){
251  //std::cout << "bas id: " << itDSBaseline->id << " raw id: " << detId << std::endl;
252  if(itDSBaseline->id != detId){
253  std::cout << "Collections out of Synch. Something of fishy is going on ;-)" << std::endl;
254  return;
255  }
256  }
257 
258 
259  actualModule_++;
260  edm::RunNumber_t const run = e.id().run();
261  edm::EventNumber_t const event = e.id().event();
262  //std::cout << "processing module N: " << actualModule_<< " detId: " << detId << " event: "<< event << std::endl;
263 
264 
265 
266  edm::DetSet<SiStripRawDigi>::const_iterator itRaw = itRawDigis->begin();
267  bool restAPV[6] = {0,0,0,0,0,0};
268  int strip =0, totADC=0;
269  int minAPVRes = 7, maxAPVRes = -1;
270  for(;itRaw != itRawDigis->end(); ++itRaw, ++strip){
271  float adc = itRaw->adc();
272  totADC+= adc;
273  if(strip%127 ==0){
274  //std::cout << "totADC " << totADC << std::endl;
275  int APV = strip/128;
276  if(totADC!= 0){
277  restAPV[APV] = true;
278  totADC =0;
279  if(APV>maxAPVRes) maxAPVRes = APV;
280  if(APV<minAPVRes) minAPVRes = APV;
281  }
282  }
283  }
284 
285  uint16_t bins =768;
286  float minx = -0.5, maxx=767.5;
287  if(minAPVRes !=7){
288  minx = minAPVRes * 128 -0.5;
289  maxx = maxAPVRes * 128 + 127.5;
290  bins = maxx-minx;
291  }
292 
293  sprintf(detIds,"%ul", detId);
294  sprintf(evs,"%llu", event);
295  sprintf(runs,"%u", run);
296  char* dHistoName = Form("Id:%s_run:%s_ev:%s",detIds, runs, evs);
297  h1ProcessedRawDigis_ = sdProcessedRawDigis_.make<TH1F>(dHistoName,dHistoName, bins, minx, maxx);
298 
299  if(plotBaseline_){
300  h1Baseline_ = sdBaseline_.make<TH1F>(dHistoName,dHistoName, bins, minx, maxx);
301  h1Baseline_->SetXTitle("strip#");
302  h1Baseline_->SetYTitle("ADC");
303  h1Baseline_->SetMaximum(1024.);
304  h1Baseline_->SetMinimum(-300.);
305  h1Baseline_->SetLineWidth(2);
306  h1Baseline_->SetLineStyle(2);
307  h1Baseline_->SetLineColor(2);
308  }
309 
310  if(plotClusters_){
311  h1Clusters_ = sdClusters_.make<TH1F>(dHistoName,dHistoName, bins, minx, maxx);
312 
313  h1Clusters_->SetXTitle("strip#");
314  h1Clusters_->SetYTitle("ADC");
315  h1Clusters_->SetMaximum(1024.);
316  h1Clusters_->SetMinimum(-300.);
317  h1Clusters_->SetLineWidth(2);
318  h1Clusters_->SetLineStyle(2);
319  h1Clusters_->SetLineColor(3);
320  }
321 
322  h1ProcessedRawDigis_->SetXTitle("strip#");
323  h1ProcessedRawDigis_->SetYTitle("ADC");
324  h1ProcessedRawDigis_->SetMaximum(1024.);
325  h1ProcessedRawDigis_->SetMinimum(-300.);
326  h1ProcessedRawDigis_->SetLineWidth(2);
327 
328  std::vector<int16_t> ProcessedRawDigis(itRawDigis->size());
329  subtractorPed_->subtract( *itRawDigis, ProcessedRawDigis);
330 
332  if(plotBaseline_) itBaseline = itDSBaseline->begin();
333  std::vector<int16_t>::const_iterator itProcessedRawDigis;
334 
335  strip =0;
336  for(itProcessedRawDigis = ProcessedRawDigis.begin();itProcessedRawDigis != ProcessedRawDigis.end(); ++itProcessedRawDigis){
337  if(restAPV[strip/128]){
338  float adc = *itProcessedRawDigis;
339  h1ProcessedRawDigis_->Fill(strip, adc);
340  if(plotBaseline_){
341  h1Baseline_->Fill(strip, itBaseline->adc());
342  ++itBaseline;
343  }
344  }
345  ++strip;
346  }
347 
348  if(plotBaseline_) ++itDSBaseline;
349  if(plotClusters_){
350  edmNew::DetSetVector<SiStripCluster>::const_iterator itClusters = clusters->begin();
351  for ( ; itClusters != clusters->end(); ++itClusters ){
352  for ( edmNew::DetSet<SiStripCluster>::const_iterator clus = itClusters->begin(); clus != itClusters->end(); ++clus){
353  if(itClusters->id() == detId){
354  int firststrip = clus->firstStrip();
355  //std::cout << "Found cluster in detId " << detId << " " << firststrip << " " << clus->amplitudes().size() << " -----------------------------------------------" << std::endl;
356  strip=0;
357  for( auto itAmpl = clus->amplitudes().begin(); itAmpl != clus->amplitudes().end(); ++itAmpl){
358  h1Clusters_->Fill(firststrip+strip, *itAmpl);
359  ++strip;
360  }
361  }
362  }
363  }
364  }
365  }
366 
367 }
368 
369 
370 // ------------ method called once each job just before starting event loop ------------
372 {
373 
374 
375 actualModule_ =0;
376 
377 
378 }
379 
380 // ------------ method called once each job just after ending the event loop ------------
381 void
383 
384 }
385 
386 //define this as a plug-in
388 
int adc(sample_type sample)
get the ADC sample (12 bits)
RunNumber_t run() const
Definition: EventID.h:39
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator end(bool update=false) const
virtual void endJob() override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
unsigned long long EventNumber_t
std::vector< TH1F > vBaselineHisto_
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
data_type const * const_iterator
Definition: DetSetNew.h:30
std::pair< ContainerIterator, ContainerIterator > Range
virtual void beginJob() override
id_type id(size_t cell) const
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
SiStripBaselineAnalyzer(const edm::ParameterSet &)
std::vector< TH1F > vBaselinePointsHisto_
T * make(const Args &...args) const
make new ROOT object
edm::Service< TFileService > fs_
std::vector< TH1F > vClusterHisto_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:416
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:361
size_type size() const
Return the number of contained DetSets.
Definition: DetSetVector.h:283
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
void allPeds(std::vector< int > &pefs, const Range &range) const
void getDetIds(std::vector< uint32_t > &DetIds_) const
const T & get() const
Definition: EventSetup.h:56
edm::EventID id() const
Definition: EventBase.h:60
HLT enums.
std::vector< TH1F > vProcessedRawDigiHisto_
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:346
unsigned int RunNumber_t
const Range getRange(const uint32_t &detID) const
edm::ESHandle< SiStripPedestals > pedestalsHandle
collection_type::const_iterator const_iterator
Definition: DetSet.h:33
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:104
std::auto_ptr< SiStripPedestalsSubtractor > subtractorPed_
const_iterator begin(bool update=false) const
static std::auto_ptr< SiStripPedestalsSubtractor > create_SubtractorPed(const edm::ParameterSet &)
Definition: event.py:1