CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripHitEffFromCalibTree.cc
Go to the documentation of this file.
1 //Original Author: Christopher Edelmaier
2 // Created: Feb. 11, 2010
3 #include <memory>
4 #include <string>
5 #include <iostream>
6 
15 
41 
54 
57 
62 
63 #include "TFile.h"
64 #include "TCanvas.h"
65 #include "TObjString.h"
66 #include "TString.h"
67 #include "TH1F.h"
68 #include "TH2F.h"
69 #include "TProfile.h"
70 #include "TF1.h"
71 #include "TROOT.h"
72 #include "TTree.h"
73 #include "TChain.h"
74 #include "TStyle.h"
75 #include "TLeaf.h"
76 #include "TGaxis.h"
77 #include "TGraphAsymmErrors.h"
78 #include "TLatex.h"
79 #include "TLegend.h"
80 
81 using namespace edm;
82 using namespace reco;
83 using namespace std;
84 
85 struct hit{
86  double x;
87  double y;
88  double z;
89  unsigned int id;
90 };
91 
92 class SiStripHitEffFromCalibTree : public ConditionDBWriter<SiStripBadStrip> {
93  public:
96 
97  private:
98  virtual void algoBeginJob();
99  virtual void algoEndJob() override;
100  virtual void algoAnalyze(const edm::Event& e, const edm::EventSetup& c) override;
101  void SetBadComponents(int i, int component,SiStripQuality::BadComponent& BC, std::stringstream ssV[4][19], int NBadComponent[4][19][4]);
102  void makeTKMap();
103  void makeHotColdMaps();
104  void makeSQLite();
105  void totalStatistics();
106  void makeSummary();
107  float calcPhi(float x, float y);
108 
113  SiStripBadStrip* getNewObject() override;
114 
116  TTree* CalibTree;
118  float threshold;
119  unsigned int nModsMin;
120  unsigned int doSummary;
121  float _ResXSig;
122  unsigned int _bunchx;
123  vector<hit> hits[23];
124  vector<TH2F*> HotColdMaps;
125  map< unsigned int, pair< unsigned int, unsigned int> > modCounter[23];
128  int layerfound[23];
129  int layertotal[23];
130  int goodlayertotal[35];
131  int goodlayerfound[35];
132  int alllayertotal[35];
133  int alllayerfound[35];
134  map< unsigned int, double > BadModules;
135 };
136 
139  FileInPath_("CalibTracker/SiStripCommon/data/SiStripDetInfo.dat")
140 {
141  CalibTreeFilename = conf.getParameter<std::string>("CalibTreeFilename");
142  threshold = conf.getParameter<double>("Threshold");
143  nModsMin = conf.getParameter<int>("nModsMin");
144  doSummary = conf.getParameter<int>("doSummary");
145  _ResXSig = conf.getUntrackedParameter<double>("ResXSig",-1);
146  _bunchx = conf.getUntrackedParameter<int>("BunchCrossing",0);
148 
149  quality_ = new SiStripQuality;
150 }
151 
153 
155  //I have no idea what goes here
156  //fs->make<TTree>("HitEffHistos","Tree of the inefficient hit histograms");
157 }
158 
160  //Still have no idea what goes here
161 
162 }
163 
165  //Retrieve tracker topology from geometry
166  edm::ESHandle<TrackerTopology> tTopoHandle;
167  c.get<IdealGeometryRecord>().get(tTopoHandle);
168  const TrackerTopology* const tTopo = tTopoHandle.product();
169 
170  //Open the ROOT Calib Tree
171  CalibTreeFile = TFile::Open(CalibTreeFilename,"READ");
172  CalibTreeFile->cd("anEff");
173  CalibTree = (TTree*)(gDirectory->Get("traj")) ;
174  TLeaf* BadLf = CalibTree->GetLeaf("ModIsBad");
175  TLeaf* sistripLf = CalibTree->GetLeaf("SiStripQualBad");
176  TLeaf* idLf = CalibTree->GetLeaf("Id");
177  TLeaf* acceptLf = CalibTree->GetLeaf("withinAcceptance");
178  TLeaf* layerLf = CalibTree->GetLeaf("layer");
179  TLeaf* nHitsLf = CalibTree->GetLeaf("nHits");
180  TLeaf* xLf = CalibTree->GetLeaf("TrajGlbX");
181  TLeaf* yLf = CalibTree->GetLeaf("TrajGlbY");
182  TLeaf* zLf = CalibTree->GetLeaf("TrajGlbZ");
183  TLeaf* ResXSigLf = CalibTree->GetLeaf("ResXSig");
184  TLeaf* BunchLf(0);
185  for(int l=0; l < 35; l++) {
186  goodlayertotal[l] = 0;
187  goodlayerfound[l] = 0;
188  alllayertotal[l] = 0;
189  alllayerfound[l] = 0;
190  }
191  if(_bunchx != 0) {
192  BunchLf = CalibTree->GetLeaf("bunchx");
193  }
194  int nevents = CalibTree->GetEntries();
195  cout << "Successfully loaded analyze function with " << nevents << " events!\n";
196  cout << "A module is bad if efficiency < " << threshold << " and has at least " << nModsMin << " nModsMin." << endl;
197 
198  //Loop through all of the events
199  for(int j =0; j < nevents; j++) {
200  CalibTree->GetEvent(j);
201  unsigned int isBad = (unsigned int)BadLf->GetValue();
202  unsigned int quality = (unsigned int)sistripLf->GetValue();
203  unsigned int id = (unsigned int)idLf->GetValue();
204  unsigned int accept = (unsigned int)acceptLf->GetValue();
205  unsigned int layer = (unsigned int)layerLf->GetValue();
206  unsigned int nHits = (unsigned int)nHitsLf->GetValue();
207  double x = xLf->GetValue();
208  double y = yLf->GetValue();
209  double z = zLf->GetValue();
210  double resxsig = ResXSigLf->GetValue();
211  bool badquality = false;
212  if(_bunchx != 0) {
213  if(_bunchx != BunchLf->GetValue()) continue;
214  }
215  //We have two things we want to do, both an XY color plot, and the efficiency measurement
216  //First, ignore anything that isn't in acceptance and isn't good quality
217 
218  //if(quality == 1 || accept != 1 || nHits < 8) continue;
219  if(accept != 1 || nHits < 8) continue;
220  if(quality == 1) badquality = true;
221 
222  //Now that we have a good event, we need to look at if we expected it or not, and the location
223  //if we didn't
224  //Fill the missing hit information first
225  bool badflag = false;
226  if(_ResXSig < 0) {
227  if(isBad == 1) badflag = true;
228  }
229  else {
230  if(isBad == 1 || resxsig > _ResXSig) badflag = true;
231  }
232  if(badflag && !badquality) {
233  hit temphit;
234  temphit.x = x;
235  temphit.y = y;
236  temphit.z = z;
237  temphit.id = id;
238  hits[layer].push_back(temphit);
239  }
240  pair<unsigned int, unsigned int> newgoodpair (1,1);
241  pair<unsigned int, unsigned int> newbadpair (1,0);
242  //First, figure out if the module already exists in the map of maps
243  map< unsigned int, pair< unsigned int, unsigned int> >::iterator it = modCounter[layer].find(id);
244  if(!badquality) {
245  if(it == modCounter[layer].end()) {
246  if(badflag) modCounter[layer][id] = newbadpair;
247  else modCounter[layer][id] = newgoodpair;
248  }
249  else {
250  ((*it).second.first)++;
251  if(!badflag) ((*it).second.second)++;
252  }
253  //Have to do the decoding for which side to go on (ugh)
254  if(layer <= 10) {
255  if(!badflag) goodlayerfound[layer]++;
256  goodlayertotal[layer]++;
257  }
258  else if(layer > 10 && layer < 14) {
259  if( ((id>>13)&0x3) == 1) {
260  if(!badflag) goodlayerfound[layer]++;
261  goodlayertotal[layer]++;
262  }
263  else if( ((id>>13)&0x3) == 2) {
264  if(!badflag) goodlayerfound[layer+3]++;
265  goodlayertotal[layer+3]++;
266  }
267  }
268  else if(layer > 13 && layer <= 22) {
269  if( ((id>>18)&0x3) == 1) {
270  if(!badflag) goodlayerfound[layer+3]++;
271  goodlayertotal[layer+3]++;
272  }
273  else if( ((id>>18)&0x3) == 2) {
274  if(!badflag) goodlayerfound[layer+12]++;
275  goodlayertotal[layer+12]++;
276  }
277  }
278  }
279  //Do the one where we don't exclude bad modules!
280  if(layer <= 10) {
281  if(!badflag) alllayerfound[layer]++;
282  alllayertotal[layer]++;
283  }
284  else if(layer > 10 && layer < 14) {
285  if( ((id>>13)&0x3) == 1) {
286  if(!badflag) alllayerfound[layer]++;
287  alllayertotal[layer]++;
288  }
289  else if( ((id>>13)&0x3) == 2) {
290  if(!badflag) alllayerfound[layer+3]++;
291  alllayertotal[layer+3]++;
292  }
293  }
294  else if(layer > 13 && layer <= 22) {
295  if( ((id>>18)&0x3) == 1) {
296  if(!badflag) alllayerfound[layer+3]++;
297  alllayertotal[layer+3]++;
298  }
299  else if( ((id>>18)&0x3) == 2) {
300  if(!badflag) alllayerfound[layer+12]++;
301  alllayertotal[layer+12]++;
302  }
303  }
304  //At this point, both of our maps are loaded with the correct information
305  }
306  //CalibTreeFile->Close();
307  makeHotColdMaps();
308  makeTKMap();
309  makeSQLite();
310  totalStatistics();
311  makeSummary();
312 
314  //try to write out what's in the quality record
316  int NTkBadComponent[4]; //k: 0=BadModule, 1=BadFiber, 2=BadApv, 3=BadStrips
317  int NBadComponent[4][19][4];
318  //legend: NBadComponent[i][j][k]= SubSystem i, layer/disk/wheel j, BadModule/Fiber/Apv k
319  // i: 0=TIB, 1=TID, 2=TOB, 3=TEC
320  // k: 0=BadModule, 1=BadFiber, 2=BadApv, 3=BadStrips
321  std::stringstream ssV[4][19];
322 
323  for(int i=0;i<4;++i){
324  NTkBadComponent[i]=0;
325  for(int j=0;j<19;++j){
326  ssV[i][j].str("");
327  for(int k=0;k<4;++k)
328  NBadComponent[i][j][k]=0;
329  }
330  }
331 
332 
333  std::vector<SiStripQuality::BadComponent> BC = quality_->getBadComponentList();
334 
335  for (size_t i=0;i<BC.size();++i){
336 
337  //&&&&&&&&&&&&&
338  //Full Tk
339  //&&&&&&&&&&&&&
340 
341  if (BC[i].BadModule)
342  NTkBadComponent[0]++;
343  if (BC[i].BadFibers)
344  NTkBadComponent[1]+= ( (BC[i].BadFibers>>2)&0x1 )+ ( (BC[i].BadFibers>>1)&0x1 ) + ( (BC[i].BadFibers)&0x1 );
345  if (BC[i].BadApvs)
346  NTkBadComponent[2]+= ( (BC[i].BadApvs>>5)&0x1 )+ ( (BC[i].BadApvs>>4)&0x1 ) + ( (BC[i].BadApvs>>3)&0x1 ) +
347  ( (BC[i].BadApvs>>2)&0x1 )+ ( (BC[i].BadApvs>>1)&0x1 ) + ( (BC[i].BadApvs)&0x1 );
348 
349  //&&&&&&&&&&&&&&&&&
350  //Single SubSystem
351  //&&&&&&&&&&&&&&&&&
352 
353  int component;
354  SiStripDetId a(BC[i].detid);
355  if ( a.subdetId() == SiStripDetId::TIB ){
356  //&&&&&&&&&&&&&&&&&
357  //TIB
358  //&&&&&&&&&&&&&&&&&
359 
360  component=tTopo->tibLayer(BC[i].detid);
361  SetBadComponents(0, component, BC[i], ssV, NBadComponent);
362 
363  } else if ( a.subdetId() == SiStripDetId::TID ) {
364  //&&&&&&&&&&&&&&&&&
365  //TID
366  //&&&&&&&&&&&&&&&&&
367 
368  component=tTopo->tidSide(BC[i].detid)==2?tTopo->tidWheel(BC[i].detid):tTopo->tidWheel(BC[i].detid)+3;
369  SetBadComponents(1, component, BC[i], ssV, NBadComponent);
370 
371  } else if ( a.subdetId() == SiStripDetId::TOB ) {
372  //&&&&&&&&&&&&&&&&&
373  //TOB
374  //&&&&&&&&&&&&&&&&&
375 
376  component=tTopo->tobLayer(BC[i].detid);
377  SetBadComponents(2, component, BC[i], ssV, NBadComponent);
378 
379  } else if ( a.subdetId() == SiStripDetId::TEC ) {
380  //&&&&&&&&&&&&&&&&&
381  //TEC
382  //&&&&&&&&&&&&&&&&&
383 
384  component=tTopo->tecSide(BC[i].detid)==2?tTopo->tecWheel(BC[i].detid):tTopo->tecWheel(BC[i].detid)+9;
385  SetBadComponents(3, component, BC[i], ssV, NBadComponent);
386 
387  }
388  }
389 
390  //&&&&&&&&&&&&&&&&&&
391  // Single Strip Info
392  //&&&&&&&&&&&&&&&&&&
393  float percentage=0;
394 
395  SiStripQuality::RegistryIterator rbegin = quality_->getRegistryVectorBegin();
396  SiStripQuality::RegistryIterator rend = quality_->getRegistryVectorEnd();
397 
398  for (SiStripBadStrip::RegistryIterator rp=rbegin; rp != rend; ++rp) {
399  unsigned int detid=rp->detid;
400 
401  int subdet=-999; int component=-999;
402  SiStripDetId a(detid);
403  if ( a.subdetId() == 3 ){
404  subdet=0;
405  component=tTopo->tibLayer(detid);
406  } else if ( a.subdetId() == 4 ) {
407  subdet=1;
408  component=tTopo->tidSide(detid)==2?tTopo->tidWheel(detid):tTopo->tidWheel(detid)+3;
409  } else if ( a.subdetId() == 5 ) {
410  subdet=2;
411  component=tTopo->tobLayer(detid);
412  } else if ( a.subdetId() == 6 ) {
413  subdet=3;
414  component=tTopo->tecSide(detid)==2?tTopo->tecWheel(detid):tTopo->tecWheel(detid)+9;
415  }
416 
417  SiStripQuality::Range sqrange = SiStripQuality::Range( quality_->getDataVectorBegin()+rp->ibegin , quality_->getDataVectorBegin()+rp->iend );
418 
419  percentage=0;
420  for(int it=0;it<sqrange.second-sqrange.first;it++){
421  unsigned int range=quality_->decode( *(sqrange.first+it) ).range;
422  NTkBadComponent[3]+=range;
423  NBadComponent[subdet][0][3]+=range;
424  NBadComponent[subdet][component][3]+=range;
425  percentage+=range;
426  }
427  if(percentage!=0)
428  percentage/=128.*reader->getNumberOfApvsAndStripLength(detid).first;
429  if(percentage>1)
430  edm::LogError("SiStripQualityStatistics") << "PROBLEM detid " << detid << " value " << percentage<< std::endl;
431  }
432  //&&&&&&&&&&&&&&&&&&
433  // printout
434  //&&&&&&&&&&&&&&&&&&
435 
436  cout << "\n-----------------\nNew IOV starting from run " << e.id().run() << " event " << e.id().event() << " lumiBlock " << e.luminosityBlock() << " time " << e.time().value() << "\n-----------------\n";
437  cout << "\n-----------------\nGlobal Info\n-----------------";
438  cout << "\nBadComponent \t Modules \tFibers \tApvs\tStrips\n----------------------------------------------------------------";
439  cout << "\nTracker:\t\t"<<NTkBadComponent[0]<<"\t"<<NTkBadComponent[1]<<"\t"<<NTkBadComponent[2]<<"\t"<<NTkBadComponent[3];
440  cout << endl;
441  cout << "\nTIB:\t\t\t"<<NBadComponent[0][0][0]<<"\t"<<NBadComponent[0][0][1]<<"\t"<<NBadComponent[0][0][2]<<"\t"<<NBadComponent[0][0][3];
442  cout << "\nTID:\t\t\t"<<NBadComponent[1][0][0]<<"\t"<<NBadComponent[1][0][1]<<"\t"<<NBadComponent[1][0][2]<<"\t"<<NBadComponent[1][0][3];
443  cout << "\nTOB:\t\t\t"<<NBadComponent[2][0][0]<<"\t"<<NBadComponent[2][0][1]<<"\t"<<NBadComponent[2][0][2]<<"\t"<<NBadComponent[2][0][3];
444  cout << "\nTEC:\t\t\t"<<NBadComponent[3][0][0]<<"\t"<<NBadComponent[3][0][1]<<"\t"<<NBadComponent[3][0][2]<<"\t"<<NBadComponent[3][0][3];
445  cout << "\n";
446 
447  for (int i=1;i<5;++i)
448  cout << "\nTIB Layer " << i << " :\t\t"<<NBadComponent[0][i][0]<<"\t"<<NBadComponent[0][i][1]<<"\t"<<NBadComponent[0][i][2]<<"\t"<<NBadComponent[0][i][3];
449  cout << "\n";
450  for (int i=1;i<4;++i)
451  cout << "\nTID+ Disk " << i << " :\t\t"<<NBadComponent[1][i][0]<<"\t"<<NBadComponent[1][i][1]<<"\t"<<NBadComponent[1][i][2]<<"\t"<<NBadComponent[1][i][3];
452  for (int i=4;i<7;++i)
453  cout << "\nTID- Disk " << i-3 << " :\t\t"<<NBadComponent[1][i][0]<<"\t"<<NBadComponent[1][i][1]<<"\t"<<NBadComponent[1][i][2]<<"\t"<<NBadComponent[1][i][3];
454  cout << "\n";
455  for (int i=1;i<7;++i)
456  cout << "\nTOB Layer " << i << " :\t\t"<<NBadComponent[2][i][0]<<"\t"<<NBadComponent[2][i][1]<<"\t"<<NBadComponent[2][i][2]<<"\t"<<NBadComponent[2][i][3];
457  cout << "\n";
458  for (int i=1;i<10;++i)
459  cout << "\nTEC+ Disk " << i << " :\t\t"<<NBadComponent[3][i][0]<<"\t"<<NBadComponent[3][i][1]<<"\t"<<NBadComponent[3][i][2]<<"\t"<<NBadComponent[3][i][3];
460  for (int i=10;i<19;++i)
461  cout << "\nTEC- Disk " << i-9 << " :\t\t"<<NBadComponent[3][i][0]<<"\t"<<NBadComponent[3][i][1]<<"\t"<<NBadComponent[3][i][2]<<"\t"<<NBadComponent[3][i][3];
462  cout << "\n";
463 
464  cout << "\n----------------------------------------------------------------\n\t\t Detid \tModules Fibers Apvs\n----------------------------------------------------------------";
465  for (int i=1;i<5;++i)
466  cout << "\nTIB Layer " << i << " :" << ssV[0][i].str();
467  cout << "\n";
468  for (int i=1;i<4;++i)
469  cout << "\nTID+ Disk " << i << " :" << ssV[1][i].str();
470  for (int i=4;i<7;++i)
471  cout << "\nTID- Disk " << i-3 << " :" << ssV[1][i].str();
472  cout << "\n";
473  for (int i=1;i<7;++i)
474  cout << "\nTOB Layer " << i << " :" << ssV[2][i].str();
475  cout << "\n";
476  for (int i=1;i<10;++i)
477  cout << "\nTEC+ Disk " << i << " :" << ssV[3][i].str();
478  for (int i=10;i<19;++i)
479  cout << "\nTEC- Disk " << i-9 << " :" << ssV[3][i].str();
480 
481 }
482 
484  cout << "Entering hot cold map generation!\n";
485  TStyle* gStyle = new TStyle("gStyle","myStyle");
486  gStyle->cd();
487  gStyle->SetPalette(1);
488  gStyle->SetCanvasColor(kWhite);
489  gStyle->SetOptStat(0);
490  //Here we make the hot/cold color maps that we love so very much
491  //Already have access to the data as a private variable
492  //Create all of the histograms in the TFileService
493  TH2F *temph2;
494  for(Long_t maplayer = 1; maplayer <=22; maplayer++) {
495  //Initialize all of the histograms
496  if(maplayer > 0 && maplayer <= 4) {
497  //We are in the TIB
498  temph2 = fs->make<TH2F>(Form("%s%i","TIB",(int)(maplayer)),"TIB",100,-1,361,100,-100,100);
499  temph2->GetXaxis()->SetTitle("Phi");
500  temph2->GetXaxis()->SetBinLabel(1,TString("360"));
501  temph2->GetXaxis()->SetBinLabel(50,TString("180"));
502  temph2->GetXaxis()->SetBinLabel(100,TString("0"));
503  temph2->GetYaxis()->SetTitle("Global Z");
504  temph2->SetOption("colz");
505  HotColdMaps.push_back(temph2);
506  }
507  else if(maplayer > 4 && maplayer <= 10) {
508  //We are in the TOB
509  temph2 = fs->make<TH2F>(Form("%s%i","TOB",(int)(maplayer-4)),"TOB",100,-1,361,100,-120,120);
510  temph2->GetXaxis()->SetTitle("Phi");
511  temph2->GetXaxis()->SetBinLabel(1,TString("360"));
512  temph2->GetXaxis()->SetBinLabel(50,TString("180"));
513  temph2->GetXaxis()->SetBinLabel(100,TString("0"));
514  temph2->GetYaxis()->SetTitle("Global Z");
515  temph2->SetOption("colz");
516  HotColdMaps.push_back(temph2);
517  }
518  else if(maplayer > 10 && maplayer <= 13) {
519  //We are in the TID
520  //Split by +/-
521  temph2 = fs->make<TH2F>(Form("%s%i","TID-",(int)(maplayer-10)),"TID-",100,-100,100,100,-100,100);
522  temph2->GetXaxis()->SetTitle("Global Y");
523  temph2->GetXaxis()->SetBinLabel(1,TString("+Y"));
524  temph2->GetXaxis()->SetBinLabel(50,TString("0"));
525  temph2->GetXaxis()->SetBinLabel(100,TString("-Y"));
526  temph2->GetYaxis()->SetTitle("Global X");
527  temph2->GetYaxis()->SetBinLabel(1,TString("-X"));
528  temph2->GetYaxis()->SetBinLabel(50,TString("0"));
529  temph2->GetYaxis()->SetBinLabel(100,TString("+X"));
530  temph2->SetOption("colz");
531  HotColdMaps.push_back(temph2);
532  temph2 = fs->make<TH2F>(Form("%s%i","TID+",(int)(maplayer-10)),"TID+",100,-100,100,100,-100,100);
533  temph2->GetXaxis()->SetTitle("Global Y");
534  temph2->GetXaxis()->SetBinLabel(1,TString("+Y"));
535  temph2->GetXaxis()->SetBinLabel(50,TString("0"));
536  temph2->GetXaxis()->SetBinLabel(100,TString("-Y"));
537  temph2->GetYaxis()->SetTitle("Global X");
538  temph2->GetYaxis()->SetBinLabel(1,TString("-X"));
539  temph2->GetYaxis()->SetBinLabel(50,TString("0"));
540  temph2->GetYaxis()->SetBinLabel(100,TString("+X"));
541  temph2->SetOption("colz");
542  HotColdMaps.push_back(temph2);
543  }
544  else if(maplayer > 13) {
545  //We are in the TEC
546  //Split by +/-
547  temph2 = fs->make<TH2F>(Form("%s%i","TEC-",(int)(maplayer-13)),"TEC-",100,-120,120,100,-120,120);
548  temph2->GetXaxis()->SetTitle("Global Y");
549  temph2->GetXaxis()->SetBinLabel(1,TString("+Y"));
550  temph2->GetXaxis()->SetBinLabel(50,TString("0"));
551  temph2->GetXaxis()->SetBinLabel(100,TString("-Y"));
552  temph2->GetYaxis()->SetTitle("Global X");
553  temph2->GetYaxis()->SetBinLabel(1,TString("-X"));
554  temph2->GetYaxis()->SetBinLabel(50,TString("0"));
555  temph2->GetYaxis()->SetBinLabel(100,TString("+X"));
556  temph2->SetOption("colz");
557  HotColdMaps.push_back(temph2);
558  temph2 = fs->make<TH2F>(Form("%s%i","TEC+",(int)(maplayer-13)),"TEC+",100,-120,120,100,-120,120);
559  temph2->GetXaxis()->SetTitle("Global Y");
560  temph2->GetXaxis()->SetBinLabel(1,TString("+Y"));
561  temph2->GetXaxis()->SetBinLabel(50,TString("0"));
562  temph2->GetXaxis()->SetBinLabel(100,TString("-Y"));
563  temph2->GetYaxis()->SetTitle("Global X");
564  temph2->GetYaxis()->SetBinLabel(1,TString("-X"));
565  temph2->GetYaxis()->SetBinLabel(50,TString("0"));
566  temph2->GetYaxis()->SetBinLabel(100,TString("+X"));
567  temph2->SetOption("colz");
568  HotColdMaps.push_back(temph2);
569  }
570  }
571  for(Long_t mylayer = 1; mylayer <= 22; mylayer++) {
572  //Determine what kind of plot we want to write out
573  //Loop through the entirety of each layer
574  //Create an array of the histograms
575  vector<hit>::const_iterator iter;
576  for(iter = hits[mylayer].begin(); iter != hits[mylayer].end(); iter++) {
577  //Looping over the particular layer
578  //Fill by 360-x to get the proper location to compare with TKMaps of phi
579  //Also global xy is messed up
580  if(mylayer > 0 && mylayer <= 4) {
581  //We are in the TIB
582  float phi = calcPhi(iter->x, iter->y);
583  HotColdMaps[mylayer - 1]->Fill(360.-phi,iter->z,1.);
584  }
585  else if(mylayer > 4 && mylayer <= 10) {
586  //We are in the TOB
587  float phi = calcPhi(iter->x,iter->y);
588  HotColdMaps[mylayer - 1]->Fill(360.-phi,iter->z,1.);
589  }
590  else if(mylayer > 10 && mylayer <= 13) {
591  //We are in the TID
592  //There are 2 different maps here
593  int side = (((iter->id)>>13) & 0x3);
594  if(side == 1) HotColdMaps[(mylayer - 1) + (mylayer - 11)]->Fill(-iter->y,iter->x,1.);
595  else if(side == 2) HotColdMaps[(mylayer - 1) + (mylayer - 10)]->Fill(-iter->y,iter->x,1.);
596  //if(side == 1) HotColdMaps[(mylayer - 1) + (mylayer - 11)]->Fill(iter->x,iter->y,1.);
597  //else if(side == 2) HotColdMaps[(mylayer - 1) + (mylayer - 10)]->Fill(iter->x,iter->y,1.);
598  }
599  else if(mylayer > 13) {
600  //We are in the TEC
601  //There are 2 different maps here
602  int side = (((iter->id)>>18) & 0x3);
603  if(side == 1) HotColdMaps[(mylayer + 2) + (mylayer - 14)]->Fill(-iter->y,iter->x,1.);
604  else if(side == 2) HotColdMaps[(mylayer + 2) + (mylayer - 13)]->Fill(-iter->y,iter->x,1.);
605  //if(side == 1) HotColdMaps[(mylayer + 2) + (mylayer - 14)]->Fill(iter->x,iter->y,1.);
606  //else if(side == 2) HotColdMaps[(mylayer + 2) + (mylayer - 13)]->Fill(iter->x,iter->y,1.);
607  }
608  }
609  }
610  cout << "Finished HotCold Map Generation\n";
611 }
612 
614  cout << "Entering TKMap generation!\n";
615  tkmap = new TrackerMap(" Detector Inefficiency ");
616  tkmapbad = new TrackerMap(" Inefficient Modules ");
617  for(Long_t i = 1; i <= 22; i++) {
618  layertotal[i] = 0;
619  layerfound[i] = 0;
620  //Loop over every layer, extracting the information from
621  //the map of the efficiencies
622  map<unsigned int, pair<unsigned int, unsigned int> >::const_iterator ih;
623  for( ih = modCounter[i].begin(); ih != modCounter[i].end(); ih++) {
624  //We should be in the layer in question, and looping over all of the modules in said layer
625  //Generate the list for the TKmap, and the bad module list
626  double myeff = (double)(((*ih).second).second)/(((*ih).second).first);
627  if ( ((((*ih).second).first) >= nModsMin) && (myeff < threshold) ) {
628  //We have a bad module, put it in the list!
629  BadModules[(*ih).first] = myeff;
630  tkmapbad->fillc((*ih).first,255,0,0);
631  cout << "Layer " << i << " module " << (*ih).first << " efficiency " << myeff << " " << (((*ih).second).second) << "/" << (((*ih).second).first) << endl;
632  }
633  else {
634  //Fill the bad list with empty results for every module
635  tkmapbad->fillc((*ih).first,255,255,255);
636  }
637  if((((*ih).second).first) < 100 ) {
638  cout << "Module " << (*ih).first << " layer " << i << " is under occupancy at " << (((*ih).second).first) << endl;
639  }
640  //Put any module into the TKMap
641  //Should call module ID, and then 1- efficiency for that module
642  //if((*ih).first == 369137820) {
643  // cout << "Module 369137820 has 1-eff of " << 1.-myeff << endl;
644  //cout << "Which is " << ((*ih).second).second << "/" << ((*ih).second).first << endl;
645  //}
646  tkmap->fill((*ih).first,1.-myeff);
647  //Find the total number of hits in the module
648  layertotal[i] += int(((*ih).second).first);
649  layerfound[i] += int(((*ih).second).second);
650  }
651  }
652  tkmap->save(true, 0, 0, "SiStripHitEffTKMap.png");
653  tkmapbad->save(true, 0, 0, "SiStripHitEffTKMapBad.png");
654  cout << "Finished TKMap Generation\n";
655 }
656 
658  //Generate the SQLite file for use in the Database of the bad modules!
659  cout << "Entering SQLite file generation!\n";
660  std::vector<unsigned int> BadStripList;
661  unsigned short NStrips;
662  unsigned int id1;
663  SiStripQuality* pQuality = new SiStripQuality;
664  //This is the list of the bad strips, use to mask out entire APVs
665  //Now simply go through the bad hit list and mask out things that
666  //are bad!
667  map< unsigned int, double >::const_iterator it;
668  for(it = BadModules.begin(); it != BadModules.end(); it++) {
669  //We need to figure out how many strips are in this particular module
670  //To Mask correctly!
671  NStrips=reader->getNumberOfApvsAndStripLength((*it).first).first*128;
672  cout << "Number of strips module " << (*it).first << " is " << NStrips << endl;
673  BadStripList.push_back(pQuality->encode(0,NStrips,0));
674  //Now compact into a single bad module
675  id1=(unsigned int)(*it).first;
676  cout << "ID1 shoudl match list of modules above " << id1 << endl;
677  quality_->compact(id1,BadStripList);
678  SiStripQuality::Range range(BadStripList.begin(),BadStripList.end());
679  quality_->put(id1,range);
680  BadStripList.clear();
681  }
682  //Fill all the bad components now
683  quality_->fillBadComponents();
684 }
685 
687  //Calculate the statistics by layer
688  int totalfound = 0;
689  int totaltotal = 0;
690  double layereff;
691  for(Long_t i=1; i<=22; i++) {
692  layereff = double(layerfound[i])/double(layertotal[i]);
693  cout << "Layer " << i << " has total efficiency " << layereff << " " << layerfound[i] << "/" << layertotal[i] << endl;
694  totalfound += layerfound[i];
695  totaltotal += layertotal[i];
696  }
697  cout << "The total efficiency is " << double(totalfound)/double(totaltotal) << endl;
698 }
699 
701  //setTDRStyle();
702 
703  int nLayers = 34;
704 
705  TH1F *found = fs->make<TH1F>("found","found",nLayers+1,0,nLayers+1);
706  TH1F *all = fs->make<TH1F>("all","all",nLayers+1,0,nLayers+1);
707  TH1F *found2 = fs->make<TH1F>("found2","found2",nLayers+1,0,nLayers+1);
708  TH1F *all2 = fs->make<TH1F>("all2","all2",nLayers+1,0,nLayers+1);
709  // first bin only to keep real data off the y axis so set to -1
710  found->SetBinContent(0,-1);
711  all->SetBinContent(0,1);
712 
713  TCanvas *c7 =new TCanvas("c7"," test ",10,10,800,600);
714  c7->SetFillColor(0);
715  c7->SetGrid();
716 
717  for (Long_t i=1; i< nLayers+1; ++i) {
718  if (i==10) i++;
719  if (i==25) i++;
720  if (i==34) break;
721 
722  cout << "Fill only good modules layer " << i << ": S = " << goodlayerfound[i] << " B = " << goodlayertotal[i] << endl;
723  if (goodlayertotal[i] > 5) {
724  found->SetBinContent(i,goodlayerfound[i]);
725  all->SetBinContent(i,goodlayertotal[i]);
726  } else {
727  found->SetBinContent(i,0);
728  all->SetBinContent(i,10);
729  }
730 
731  cout << "Filling all modules layer " << i << ": S = " << alllayerfound[i] << " B = " << alllayertotal[i] << endl;
732  if (alllayertotal[i] > 5) {
733  found2->SetBinContent(i,alllayerfound[i]);
734  all2->SetBinContent(i,alllayertotal[i]);
735  } else {
736  found2->SetBinContent(i,0);
737  all2->SetBinContent(i,10);
738  }
739 
740  }
741 
742  found->Sumw2();
743  all->Sumw2();
744 
745  found2->Sumw2();
746  all2->Sumw2();
747 
748  TGraphAsymmErrors *gr = new TGraphAsymmErrors(nLayers+1);
749  gr->BayesDivide(found,all);
750 
751  TGraphAsymmErrors *gr2 = new TGraphAsymmErrors(nLayers+1);
752  gr2->BayesDivide(found2,all2);
753 
754  for(int j = 0; j<nLayers+1; j++){
755  gr->SetPointError(j, 0., 0., gr->GetErrorYlow(j),gr->GetErrorYhigh(j) );
756  gr2->SetPointError(j, 0., 0., gr2->GetErrorYlow(j),gr2->GetErrorYhigh(j) );
757  }
758 
759  gr->GetXaxis()->SetLimits(0,nLayers);
760  gr->SetMarkerColor(2);
761  gr->SetMarkerSize(1.2);
762  gr->SetLineColor(2);
763  gr->SetLineWidth(4);
764  gr->SetMarkerStyle(20);
765  gr->SetMinimum(0.90);
766  gr->SetMaximum(1.001);
767  gr->GetYaxis()->SetTitle("Efficiency");
768 
769  gr2->GetXaxis()->SetLimits(0,nLayers);
770  gr2->SetMarkerColor(1);
771  gr2->SetMarkerSize(1.2);
772  gr2->SetLineColor(1);
773  gr2->SetLineWidth(4);
774  gr2->SetMarkerStyle(21);
775  gr2->SetMinimum(0.90);
776  gr2->SetMaximum(1.001);
777  gr2->GetYaxis()->SetTitle("Efficiency");
778  //cout << "starting labels" << endl;
779  //for ( int k=1; k<nLayers+1; k++) {
780  for ( Long_t k=1; k<nLayers+1; k++) {
781  if (k==10) k++;
782  if (k==25) k++;
783  if (k==34) break;
784  TString label;
785  if (k<5) {
786  label = TString("TIB ") + k;
787  } else if (k>4&&k<11) {
788  label = TString("TOB ")+(k-4);
789  } else if (k>10&&k<14) {
790  label = TString("TID- ")+(k-10);
791  } else if (k>13&&k<17) {
792  label = TString("TID+ ")+(k-13);
793  } else if (k>16&&k<26) {
794  label = TString("TEC- ")+(k-16);
795  } else if (k>25) {
796  label = TString("TEC+ ")+(k-25);
797  }
798  gr->GetXaxis()->SetBinLabel(((k+1)*100)/(nLayers)-2,label);
799  gr2->GetXaxis()->SetBinLabel(((k+1)*100)/(nLayers)-2,label);
800  }
801 
802  gr->Draw("AP");
803  gr->GetXaxis()->SetNdivisions(36);
804 
805  c7->cd();
806  TPad *overlay = new TPad("overlay","",0,0,1,1);
807  overlay->SetFillStyle(4000);
808  overlay->SetFillColor(0);
809  overlay->SetFrameFillStyle(4000);
810  overlay->Draw("same");
811  overlay->cd();
812  gr2->Draw("AP");
813 
814  TLegend *leg = new TLegend(0.70,0.20,0.92,0.39);
815  leg->AddEntry(gr,"Good Modules","p");
816  leg->AddEntry(gr2,"All Modules","p");
817  leg->SetTextSize(0.020);
818  leg->SetFillColor(0);
819  leg->Draw("same");
820 
821  c7->SaveAs("Summary.png");
822 }
823 
825  //Need this for a Condition DB Writer
826  //Initialize a return variable
828 
829  SiStripBadStrip::RegistryIterator rIter=quality_->getRegistryVectorBegin();
830  SiStripBadStrip::RegistryIterator rIterEnd=quality_->getRegistryVectorEnd();
831 
832  for(;rIter!=rIterEnd;++rIter){
833  SiStripBadStrip::Range range(quality_->getDataVectorBegin()+rIter->ibegin,quality_->getDataVectorBegin()+rIter->iend);
834  if ( ! obj->put(rIter->detid,range) )
835  edm::LogError("SiStripHitEffFromCalibTree")<<"[SiStripHitEffFromCalibTree::getNewObject] detid already exists"<<std::endl;
836  }
837 
838  return obj;
839 }
840 
842  float phi = 0;
843  float Pi = 3.14159;
844  if((x>=0)&&(y>=0)) phi = atan(y/x);
845  else if((x>=0)&&(y<=0)) phi = atan(y/x) + 2*Pi;
846  else if((x<=0)&&(y>=0)) phi = atan(y/x) + Pi;
847  else phi = atan(y/x) + Pi;
848  phi = phi*180.0/Pi;
849 
850  return phi;
851 }
852 
853 void SiStripHitEffFromCalibTree::SetBadComponents(int i, int component,SiStripQuality::BadComponent& BC, std::stringstream ssV[4][19], int NBadComponent[4][19][4]){
854 
855  int napv=reader->getNumberOfApvsAndStripLength(BC.detid).first;
856 
857  ssV[i][component] << "\n\t\t "
858  << BC.detid
859  << " \t " << BC.BadModule << " \t "
860  << ( (BC.BadFibers)&0x1 ) << " ";
861  if (napv==4)
862  ssV[i][component] << "x " <<( (BC.BadFibers>>1)&0x1 );
863 
864  if (napv==6)
865  ssV[i][component] << ( (BC.BadFibers>>1)&0x1 ) << " "
866  << ( (BC.BadFibers>>2)&0x1 );
867  ssV[i][component] << " \t "
868  << ( (BC.BadApvs)&0x1 ) << " "
869  << ( (BC.BadApvs>>1)&0x1 ) << " ";
870  if (napv==4)
871  ssV[i][component] << "x x " << ( (BC.BadApvs>>2)&0x1 ) << " "
872  << ( (BC.BadApvs>>3)&0x1 );
873  if (napv==6)
874  ssV[i][component] << ( (BC.BadApvs>>2)&0x1 ) << " "
875  << ( (BC.BadApvs>>3)&0x1 ) << " "
876  << ( (BC.BadApvs>>4)&0x1 ) << " "
877  << ( (BC.BadApvs>>5)&0x1 ) << " ";
878 
879  if (BC.BadApvs){
880  NBadComponent[i][0][2]+= ( (BC.BadApvs>>5)&0x1 )+ ( (BC.BadApvs>>4)&0x1 ) + ( (BC.BadApvs>>3)&0x1 ) +
881  ( (BC.BadApvs>>2)&0x1 )+ ( (BC.BadApvs>>1)&0x1 ) + ( (BC.BadApvs)&0x1 );
882  NBadComponent[i][component][2]+= ( (BC.BadApvs>>5)&0x1 )+ ( (BC.BadApvs>>4)&0x1 ) + ( (BC.BadApvs>>3)&0x1 ) +
883  ( (BC.BadApvs>>2)&0x1 )+ ( (BC.BadApvs>>1)&0x1 ) + ( (BC.BadApvs)&0x1 );
884  }
885  if (BC.BadFibers){
886  NBadComponent[i][0][1]+= ( (BC.BadFibers>>2)&0x1 )+ ( (BC.BadFibers>>1)&0x1 ) + ( (BC.BadFibers)&0x1 );
887  NBadComponent[i][component][1]+= ( (BC.BadFibers>>2)&0x1 )+ ( (BC.BadFibers>>1)&0x1 ) + ( (BC.BadFibers)&0x1 );
888  }
889  if (BC.BadModule){
890  NBadComponent[i][0][0]++;
891  NBadComponent[i][component][0]++;
892  }
893 }
894 
RunNumber_t run() const
Definition: EventID.h:42
const double Pi
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
unsigned int tibLayer(const DetId &id) const
const std::pair< unsigned short, double > getNumberOfApvsAndStripLength(uint32_t detId) const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:59
unsigned int tidWheel(const DetId &id) const
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
void SetBadComponents(int i, int component, SiStripQuality::BadComponent &BC, std::stringstream ssV[4][19], int NBadComponent[4][19][4])
Registry::const_iterator RegistryIterator
SiStripDetInfoFileReader * reader
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:26
float float float z
U second(std::pair< T, U > const &p)
edm::Service< TFileService > fs
unsigned int tidSide(const DetId &id) const
void save(bool print_total=true, float minval=0., float maxval=0., std::string s="svgmap.svg", int width=1500, int height=800)
Definition: TrackerMap.cc:698
int j
Definition: DBlmapReader.cc:9
int nevents
#define end
Definition: vmac.h:37
bool first
Definition: L1TdeRCT.cc:75
SiStripBadStrip * getNewObject() override
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
map< unsigned int, pair< unsigned int, unsigned int > > modCounter[23]
tuple conf
Definition: dbtoconf.py:185
void fillc(int idmod, int RGBcode)
Definition: TrackerMap.h:103
int k[5][pyjets_maxn]
unsigned int id
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
PixelRecoRange< float > Range
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
edm::EventID id() const
Definition: EventBase.h:56
#define begin
Definition: vmac.h:30
virtual void algoAnalyze(const edm::Event &e, const edm::EventSetup &c) override
double a
Definition: hdecay.h:121
std::pair< ContainerIterator, ContainerIterator > Range
tuple cout
Definition: gather_cfg.py:121
map< unsigned int, double > BadModules
std::string fullPath() const
Definition: FileInPath.cc:165
Definition: DDAxes.h:10
bool put(const uint32_t &detID, const InputVector &vect)
unsigned int tecWheel(const DetId &id) const
TimeValue_t value() const
Definition: Timestamp.h:56
SiStripHitEffFromCalibTree(const edm::ParameterSet &)
edm::Timestamp time() const
Definition: EventBase.h:57
void fill(int layer, int ring, int nmod, float x)
Definition: TrackerMap.cc:2776
unsigned int tobLayer(const DetId &id) const
unsigned int tecSide(const DetId &id) const
Definition: DDAxes.h:10