AliPhysics  ec707b8 (ec707b8)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups Pages
BadChannelAnalysis.C
Go to the documentation of this file.
1 
20 #if !defined(__CINT__) || defined(__MAKECINT__)
21 #include <TFile.h>
22 #include <TH1F.h>
23 #include <TH2F.h>
24 #include <TH3D.h>
25 #include <TLine.h>
26 #include <Riostream.h>
27 #include <TCanvas.h>
28 #include <TGraphErrors.h>
29 #include <TGrid.h>
30 #include <TStyle.h>
31 #include <TFileMerger.h>
32 #include <TMultiGraph.h>
33 #include <TROOT.h>
34 #include <TLegend.h>
35 #include <TString.h>
36 #include <TGridCollection.h>
37 #include <TGridResult.h>
38 #include <TClonesArray.h>
39 #include <TObjString.h>
40 #include <stdio.h>
41 #include <fstream>
42 #include <iostream>
43 #endif
44 using namespace std;
45 
46 //_________________________________________________________________________
47 //_________________________________________________________________________
48 
53 void Draw(Int_t cell[], Int_t iBC, Int_t nBC, const Int_t cellref=151)
54 {
55  gROOT ->SetStyle("Plain");
56  gStyle->SetOptStat(0);
57  gStyle->SetFillColor(kWhite);
58  gStyle->SetTitleFillColor(kWhite);
59  gStyle->SetPalette(1);
60 
61  char out[120]; char title[100]; char name[100];
62  TString slide = "Cells ";
63  slide+=cell[iBC];
64  slide+="-";
65  slide+=cell[iBC+nBC-1];
66  sprintf(out,"%d-%d.gif",cell[iBC],cell[iBC+nBC-1]);
67 
68  TH2 *hCellAmplitude = (TH2*) gFile->Get("hCellAmplitude");
69  TH1 *hCellref = hCellAmplitude->ProjectionY("badcells",cellref+1,cellref+1);
70 
71  Int_t i;
72  TCanvas *c1 = new TCanvas("badcells","badcells",1000,750) ;
73  c1->SetLogy();
74  if (nBC > 6) c1->Divide(3,3) ;
75  else if (nBC > 3) c1->Divide(3,2) ;
76  else c1->Divide(3,1);
77  // hCellref->Rebin(3);
78 
79  TLegend *leg = new TLegend(0.7, 0.7, 0.9, 0.9);
80  for(i=0; i<nBC ; i++){
81  sprintf(name,"Cell %d",cell[iBC+i]) ;
82  TH1 *hCell = hCellAmplitude->ProjectionY(name,cell[iBC+i]+1,cell[iBC+i]+1);
83  c1->cd(i%9 + 1) ;
84  c1->cd(i%9 + 1)->SetLogy();
85  sprintf(title,"Cell %d Entries : %d",cell[iBC+i], (Int_t)hCell->GetEntries()) ;
86  hCell->SetLineColor(2) ;
87  // cout<<title<<endl ;
88  hCell->SetMaximum(1e5);
89  // hCell->Rebin(3);
90  hCell->SetAxisRange(0.,4.);
91  hCell->GetXaxis()->SetTitle("E (GeV)");
92  hCell->GetYaxis()->SetTitle("N Entries");
93  hCellref->SetAxisRange(0.,4.);
94  hCell->SetLineWidth(1) ;
95  hCellref->SetLineWidth(1) ;
96  hCell->SetTitle(title);
97  hCellref->SetLineColor(1) ;
98  if(i==0){
99  leg->AddEntry(hCellref,"reference","l");
100  leg->AddEntry(hCell,"current","l");;
101  }
102  hCell->Draw() ;
103  hCellref->Draw("same") ;
104  leg->Draw();
105  }
106 
107  //CREATE A PDF FILE
108  TString PdfName = "BadCellsCandidate.pdf";
109  if(nBC<9) {
110  PdfName +=")";
111  c1->Print(PdfName.Data());}
112  else if(iBC==0) {
113  PdfName +="(";
114  c1->Print(PdfName.Data());}
115  else c1->Print(PdfName.Data());
116 
117 
118  // c1->SaveAs(out);
119  delete hCellref ;
120  delete c1 ;
121  delete leg ;
122 }
123 
124 //_________________________________________________________________________
125 //_________________________________________________________________________
126 
132 void Convert(TString fCalorimeter = "EMCAL", TString period = "LHC11h", TString pass = "pass1_HLT", TString trigger= "default")
133 {
134  TH2F * hCellAmplitude = new TH2F("hCellAmplitude","Cell Amplitude",11520,0,11520,200,0,10);
135  TH1D * hNEventsProcessedPerRun = new TH1D("hNEventsProcessedPerRun","Number of processed events vs run number",200000,100000,300000);
136  FILE * pFile;
137  TString file = "/scratch/alicehp2/mas/analyse/QA/"+period+"/"+ pass + "/runlistMB.txt" ;
138  cout<<file<<endl;
139  pFile = fopen(file.Data(), "r"); //open the text file where include the run list and correct run index
140 
141  cout<<file<<endl;
142  cout << " fcalo: " << fCalorimeter << "; period: " << period << "; pass: " << pass << " trigger "<<trigger<< endl;
143 
144  Int_t ix,iy;
145  Int_t Nentr;
146  Int_t p;
147  Int_t q;
148  Int_t ncols;
149  Int_t nlines = 0 ;
150  Int_t RunId[500] ;
151 
152  Double_t x[500] ;
153  Double_t xrun[500] ;
154 
155  while (1){
156  ncols = fscanf(pFile,"%d %d ",&p,&q);
157  if (ncols< 0) break;
158  x[nlines]=p;
159  RunId[nlines]=q;
160  xrun[nlines]=1.*q;
161  nlines++;
162  }
163 
164  fclose(pFile);
165 
166  const Int_t nRun = nlines ;
167  Double_t content;
168  TString base ;
169  TString BCfile ;
170  TString direct = "CaloQA_";
171  direct += trigger;
172 
173  for(Int_t i = 0 ; i < nRun ; i++) {
174  base = "/scratch/alicehp2/germain/QA/";
175  BCfile = base + period ;
176  BCfile += trigger ;
177  BCfile += ".root";
178  base += period ;
179  base += "/";
180  base += pass ;
181  base += "/";
182  base += RunId[i] ;
183  TString infile ;
184  infile = base + ".root" ;
185  TFile *f = TFile::Open(infile);
186  base += "/" ;
187  base += trigger ;
188  TDirectoryFile *dir = (TDirectoryFile *)f->Get(direct);
189  TList *outputList = (TList*)dir->Get(direct);
190  TH2F *hAmpId;
191  TH2F *hNEvents;
192  hAmpId =(TH2F *)outputList->FindObject("EMCAL_hAmpId");
193  hNEvents =(TH2F *)outputList->FindObject("hNEvents");
194  Nentr = (Int_t)hNEvents->GetEntries();
195  if (Nentr<100) continue ;
196  hNEventsProcessedPerRun->SetBinContent(RunId[i]-100000,(Double_t)Nentr);
197 
198  for(ix=1;ix<=200;ix++){
199  for(iy=1; iy<=11520; iy++){
200  content = 0.0 ;
201  content += hAmpId->GetBinContent(ix,iy);
202  content += hCellAmplitude->GetBinContent(iy,ix);
203  if(content > 0.5)
204  hCellAmplitude->SetBinContent(iy,ix,content);
205  }
206  }
207 
208  //cout<<i<<endl;
209  if(i==0){ cout<<"Merging/Converting procedure ..." ; cout.flush();}
210  else { cout<<"..." ; cout.flush();}
211  outputList->Delete();
212  dir->Delete();
213  f->Close();
214  delete f;
215  }
216 
217 
218  TFile *BCF = TFile::Open(BCfile,"recreate");
219  hNEventsProcessedPerRun->Write();
220  hCellAmplitude->Write();
221  BCF->Close();
222  cout<<"DONE !"<<endl;
223 }
224 
225 //_________________________________________________________________________
226 //_________________________________________________________________________
227 
237 void Process(Int_t *pflag[11520][7], TH1* inhisto, Double_t Nsigma = 4.,
238  Int_t dnbins = 200, Double_t dmaxval = -1.)
239 {
240  Int_t crit = *pflag[0][0] ; //identify the criterum processed
241  Double_t goodmax= 0. ;
242  Double_t goodmin= 0. ;
243  *pflag[0][0] =1;
244 
245  if (dmaxval < 0.)
246  {
247  dmaxval = inhisto->GetMaximum()*1.01; // 1.01 - to see the last bin
248  if(crit==2 && dmaxval > 1) dmaxval =1. ;
249  }
250 
251  TH1 *distrib = new TH1F(Form("%sDistr",inhisto->GetName()), "", dnbins, inhisto->GetMinimum(), dmaxval);
252  distrib->SetXTitle(inhisto->GetYaxis()->GetTitle());
253  distrib->SetYTitle("Entries");
254 
255  // fill distribution
256  for (Int_t c = 1; c <= inhisto->GetNbinsX(); c++)
257  distrib->Fill(inhisto->GetBinContent(c));
258 
259  // draw histogram + distribution
260  TCanvas *c1 = new TCanvas(inhisto->GetName(),inhisto->GetName(), 800,400);
261  c1->Divide(2,1);
262 
263  c1->cd(1);
264  gPad->SetLeftMargin(0.14);
265  gPad->SetRightMargin(0.06);
266  gPad->SetLogy();
267  inhisto->SetTitleOffset(1.7,"Y");
268  inhisto->Draw();
269 
270  c1->cd(2);
271  gPad->SetLeftMargin(0.14);
272  gPad->SetRightMargin(0.10);
273  gPad->SetLogy();
274  distrib->Draw();
275 
276  Int_t higherbin=0,i;
277  for (i = 2; i <= distrib->GetNbinsX(); i++){
278  if(distrib->GetBinContent(higherbin) < distrib->GetBinContent(i))
279  higherbin = i ;
280  }
281 
282  for(i=higherbin ; i<=dnbins ; i++)
283  if(distrib->GetBinContent(i)<2) break ;
284  goodmax = distrib->GetBinCenter(i);
285 
286  for(i=higherbin ; i>0 ; i--)
287  if(distrib->GetBinContent(i)<2) break ;
288  goodmin = distrib->GetBinLowEdge(i);
289 
290  TF1 *fit2 = new TF1("fit2", "gaus");
291 
292  distrib->Fit(fit2, "0LQEM", "", goodmin, goodmax);
293  Double_t sig, mean, chi2ndf;
294  mean = fit2->GetParameter(1);
295  sig = fit2->GetParameter(2);
296  chi2ndf = fit2->GetChisquare()/fit2->GetNDF();
297  goodmin = mean - Nsigma*sig ;
298  goodmax = mean + Nsigma*sig ;
299 
300  // lines
301  TLine *lline = new TLine(goodmin, 0, goodmin, distrib->GetMaximum());
302  lline->SetLineColor(kOrange);
303  lline->SetLineStyle(7);
304  lline->Draw();
305 
306  TLine *rline = new TLine(goodmax, 0, goodmax, distrib->GetMaximum());
307  rline->SetLineColor(kOrange);
308  rline->SetLineStyle(7);
309  rline->Draw();
310 
311  // legend
312  TLegend *leg = new TLegend(0.60,0.82,0.9,0.88);
313  leg->AddEntry(lline, "Good region boundary","l");
314  leg->Draw("same");
315  fit2->Draw("same");
316 
317  c1->Update();
318  TString name = "criteria-" ;
319  name+= crit;
320  name+= ".gif";
321 
322  c1->SaveAs(name);
323 
324  Int_t ntot = 0, cel;
325 
326  for (Int_t c = 1; c <= inhisto->GetNbinsX(); c++) {
327  cel=(Int_t)(inhisto->GetBinLowEdge(c)+0.1);
328  if (inhisto->GetBinContent(c) < goodmin) {
329  ntot++;
330  *pflag[cel][crit]=0;
331  }
332  else if (inhisto->GetBinContent(c) > goodmax) {
333  ntot++;
334  *pflag[cel][crit]=2;
335  }
336  }
337 }
338 
339 //_________________________________________________________________________
340 //_________________________________________________________________________
341 
350 void TestCellEandN(Int_t *pflag[11520][7], Double_t Emin = 0.1, Double_t Nsigma = 4., char* hname = "hCellAmplitude", Int_t dnbins = 200)
351 {
352  TH2 *hCellAmplitude = (TH2*) gFile->Get(hname);
353 
354  // binning parameters
355  Int_t ncells = hCellAmplitude->GetNbinsX();
356  Double_t amin = hCellAmplitude->GetXaxis()->GetXmin();
357  Double_t amax = hCellAmplitude->GetXaxis()->GetXmax();
358 
359  TH1* hCellEtotal = new TH1F(Form("%s_hCellEtotal_E%.2f",hname,Emin),
360  Form("Total deposited energy, E > %.2f GeV",Emin), ncells,amin,amax);
361  hCellEtotal->SetXTitle("AbsId");
362  hCellEtotal->SetYTitle("Energy, GeV");
363 
364  TH1F *hCellNtotal = new TH1F(Form("%s_hCellNtotal_E%.2f",hname,Emin),
365  Form("Number of entries per events, E > %.2f GeV",Emin), ncells,amin,amax);
366  hCellNtotal->SetXTitle("AbsId");
367  hCellNtotal->SetYTitle("Entries");
368 
369  TH1F *hCellEtoNtotal = new TH1F(Form("%s_hCellEtoNtotal_E%.2f",hname,Emin),
370  Form("Average energy per hit, E > %.2f GeV",Emin), ncells,amin,amax);
371  hCellEtoNtotal->SetXTitle("AbsId");
372  hCellEtoNtotal->SetYTitle("Energy, GeV");
373 
374  TH1* hNEventsProcessedPerRun = (TH1*) gFile->Get("hNEventsProcessedPerRun");
375  Double_t totalevents = hNEventsProcessedPerRun->Integral(1, hNEventsProcessedPerRun->GetNbinsX());
376 
377  // fill cells
378  for (Int_t c = 1; c <= ncells; c++) {
379  Double_t Esum = 0;
380  Double_t Nsum = 0;
381 
382 
383  for (Int_t j = 1; j <= hCellAmplitude->GetNbinsY(); j++) {
384  Double_t E = hCellAmplitude->GetYaxis()->GetBinCenter(j);
385  Double_t N = hCellAmplitude->GetBinContent(c, j);
386  if (E < Emin) continue;
387  Esum += E*N;
388  Nsum += N;
389  }
390 
391  hCellEtotal->SetBinContent(c, Esum);
392  hCellNtotal->SetBinContent(c, Nsum/totalevents);
393 
394  if (Nsum > 0.5) // number of entries >= 1
395  hCellEtoNtotal->SetBinContent(c, Esum/Nsum);
396 
397  }
398 
399  delete hCellAmplitude;
400 
401  // Process(hCellEtotal, dnbins );
402  if(*pflag[0][0]==1)
403  Process(pflag, hCellEtoNtotal, Nsigma, dnbins );
404  if(*pflag[0][0]==2)
405  Process(pflag, hCellNtotal, Nsigma, dnbins );
406 }
407 
408 //_________________________________________________________________________
409 //_________________________________________________________________________
410 
415 void TestCellShapes(Int_t *pflag[11520][7], Double_t fitEmin, Double_t fitEmax, Double_t Nsigma =4., char* hname = "hCellAmplitude", Int_t dnbins = 1000)
416 {
417  TH2 *hCellAmplitude = (TH2*) gFile->Get(Form("%s",hname));
418 
419  // binning parameters
420  Int_t ncells = hCellAmplitude->GetNbinsX();
421  Double_t amin = hCellAmplitude->GetXaxis()->GetXmin();
422  Double_t amax = hCellAmplitude->GetXaxis()->GetXmax();
423 
424  // initialize histograms
425  TH1 *hFitA = new TH1F(Form("hFitA_%s",hname),"Fit A value", ncells,amin,amax);
426  hFitA->SetXTitle("AbsId");
427  hFitA->SetYTitle("A");
428 
429  TH1 *hFitB = new TH1F(Form("hFitB_%s",hname),"Fit B value", ncells,amin,amax);
430  hFitB->SetXTitle("AbsId");
431  hFitB->SetYTitle("B");
432 
433  TH1 *hFitChi2Ndf = new TH1F(Form("hFitChi2Ndf_%s",hname),"Fit #chi^{2}/ndf value", ncells,amin,amax);
434  hFitChi2Ndf->SetXTitle("AbsId");
435  hFitChi2Ndf->SetYTitle("#chi^{2}/ndf");
436 
437  Double_t maxval1=0., maxval2=0., maxval3=0.;
438  Double_t prev=0., MSA=0., AvA = 0. ; //those param are used to automaticaly determined a reasonable maxval1
439  Double_t prev2=0., MSB=0., AvB = 0. ; //those param are used to automaticaly determined a reasonable maxval2
440  Double_t prev3=0., MSki2=0., Avki2 = 0. ; //those param are used to automaticaly determined a reasonable maxval3
441  Double_t ki2=0.0 ;
442 
443  for (Int_t k = 1; k <= ncells; k++) {
444  TF1 *fit = new TF1("fit", "[0]*exp(-[1]*x)/x^2");
445  TH1 *hCell = hCellAmplitude->ProjectionY("",k,k);
446  if (hCell->GetEntries() == 0) continue;
447  // hCell->Rebin(3);
448  hCell->Fit(fit, "0QEM", "", fitEmin, fitEmax);
449  delete hCell;
450 
451  if(fit->GetParameter(0) < 5000.){
452  hFitA->SetBinContent(k, fit->GetParameter(0));
453  if(k<3000) {
454  AvA += fit->GetParameter(0);
455  if(k==2999) maxval1 = AvA/3000. ;
456  if (prev < fit->GetParameter(0)) MSA += fit->GetParameter(0) - prev;
457  else MSA -= (fit->GetParameter(0) - prev) ;
458  prev = fit->GetParameter(0);
459  }
460 
461  else
462  {
463 
464  if((fit->GetParameter(0) - maxval1) > 0. && (fit->GetParameter(0) - maxval1) < (MSA/1000.))
465  {
466  maxval1 = fit->GetParameter(0);
467  }
468  }
469  }
470  else hFitA->SetBinContent(k, 5000.);
471 
472 
473 
474  if(fit->GetParameter(1) < 5000.){
475  hFitB->SetBinContent(k, fit->GetParameter(1));
476  if(k<3000) {
477  AvB += fit->GetParameter(1);
478  if(k==2999) maxval2 = AvB/3000. ;
479  if (prev2 < fit->GetParameter(1)) MSB += fit->GetParameter(1) - prev2;
480  else MSB -= (fit->GetParameter(1) - prev2) ;
481  prev2 = fit->GetParameter(1);
482  }
483 
484  else
485  {
486 
487  if((fit->GetParameter(1) - maxval2) > 0. && (fit->GetParameter(1) - maxval2) < (MSB/1000.))
488  {
489  maxval2 = fit->GetParameter(1);
490  }
491  }
492  }
493  else hFitB->SetBinContent(k, 5000.);
494 
495 
496  if (fit->GetNDF() != 0 ) ki2 = fit->GetChisquare()/fit->GetNDF();
497  else ki2 = 1000.;
498 
499  if(ki2 < 1000.){
500  hFitChi2Ndf->SetBinContent(k, ki2);
501  if(k<3000) {
502  Avki2 += ki2;
503  if(k==2999) maxval3 = Avki2/3000. ;
504  if (prev3 < ki2) MSki2 += ki2 - prev3;
505  else MSki2 -= (ki2 - prev3) ;
506  prev3 = ki2;
507  }
508 
509  else
510  {
511 
512  if((ki2 - maxval3) > 0. && (ki2 - maxval3) < (MSki2/1000.))
513  {
514  maxval3 = ki2;
515  }
516  }
517  }
518  else hFitChi2Ndf->SetBinContent(k, 1000.);
519 
520 
521  delete fit ;
522  }
523 
524  delete hCellAmplitude;
525 
526  // if you have problem with automatic parameter :
527  // maxval1 =
528  // maxval2 =
529  // maxval3 =
530 
531 
532  if(*pflag[0][0]==3)
533  Process(pflag, hFitChi2Ndf, Nsigma, dnbins, maxval3);
534 
535 
536  if(*pflag[0][0]==4)
537  Process(pflag, hFitA, Nsigma, dnbins, maxval1);
538 
539 
540  if(*pflag[0][0]==5)
541  Process(pflag, hFitB, Nsigma, dnbins, maxval2);
542 }
543 
544 //_________________________________________________________________________
545 //_________________________________________________________________________
546 
550 void ExcludeCells(Int_t *pexclu[11520])
551 {
552  TH2 *hCellAmplitude = (TH2*) gFile->Get("hCellAmplitude");
553 
554 
555  for (Int_t c = 1; c <= 11520; c++) {
556  Double_t Nsum = 0;
557 
558  for (Int_t l = 1; l <= hCellAmplitude->GetNbinsY(); l++) {
559  Double_t N = hCellAmplitude->GetBinContent(c, l);
560  Nsum += N;
561  }
562  if(Nsum < 0.5 && *pexclu[c-1]!=1) *pexclu[c-1]=1; //trick for criterum 7
563  //if(Nsum < 0.5 ) *pexclu[c-1]=1;
564  else *pexclu[c-1]=0;
565  }
566 
567  delete hCellAmplitude;
568 }
569 
570 //_________________________________________________________________________
571 //_________________________________________________________________________
572 
576 void KillCells(Int_t filter[], Int_t nbc)
577 {
578  TH2 *hCellAmplitude = (TH2*) gFile->Get("hCellAmplitude");
579 
580  for(Int_t i =0; i<nbc; i++){
581  for(Int_t j=0; j<= hCellAmplitude->GetNbinsY() ;j++){
582  hCellAmplitude->SetBinContent(filter[i]+1,j,0) ; }}
583 
584  TH1* hNEventsProcessedPerRun = (TH1*) gFile->Get("hNEventsProcessedPerRun");
585 
586  TFile *tf = new TFile("filter.root","recreate");
587  hCellAmplitude->Write();
588  hNEventsProcessedPerRun->Write();
589  tf->Write();
590  tf->Close();
591  delete hCellAmplitude; delete hNEventsProcessedPerRun;
592 }
593 
594 //_________________________________________________________________________
595 //_________________________________________________________________________
596 
608 void PeriodAnalysis(Int_t criterum=7, Double_t Nsigma = 4.0, Double_t Emin=0.1,
609  Double_t Emax=1.0, TString file ="none")
610 {
611  Int_t newBC[11520]; // newBC[0] donne l'id de la premiere BC trouvĂ©e
612  Int_t *pexclu[11520] ;
613  Int_t exclu[11520];
614  Int_t *pflag[11520][7] ;
615  Int_t flag[11520][7];
616  Int_t bad[1000] ;
617  Int_t i, j, nb=0;
618 
619  //INIT
620  TString output, bilan;
621  if(criterum == 7) bilan = "Results.txt" ;
622  output.Form("Criterum-%d_Emin-%.2f.txt",criterum,Emin);
623  for(i=0;i<11520;i++) { exclu[i]=0; pexclu[i] =&exclu[i];
624  for(j=0;j<7;j++) { flag[i][j] =1 ; pflag[i][j] = &flag[i][j];}}
625  flag[0][0]=criterum ; //to identify the criterum tested
626 
627 
628  //CELLS EXCLUDED
629  ExcludeCells(pexclu); //exclude cells from analysis (will not appear in results)
630  if(criterum < 7){
631  cout<<"Excluded/dead cells : "<<endl;
632  for(i=0;i<11520;i++) {if(exclu[i]!=0) {cout<<i<<", " ; nb++;}}
633  cout<<"("<<nb<<")"<<endl; nb=0;}
634 
635 
636  //CRITERUM 7 : FINAL RESULT
637  if(criterum ==7) {
638  cout<<"FINAL RESULTS"<<endl;
639  ofstream fichier(bilan, ios::out | ios::trunc);
640  if(fichier){
641  fichier<<"Dead cells : "<<endl;
642  cout<<"Dead cells : "<<endl;
643  for(i=0;i<11520;i++) {
644  if(exclu[i]!=0) {fichier<<i<<", " ; cout<<i<<", " ; nb++;}}
645  fichier<<"("<<nb<<")"<<endl; cout<<"("<<nb<<")"<<endl; nb=0;
646 
647  TFile::Open("filter.root");
648  ExcludeCells(pexclu);
649  fichier<<"Bad cells candidates : "<<endl; cout<<"Bad cells candidates : "<<endl;
650  for(i=0;i<11520;i++) {
651  if(exclu[i]!=0) {bad[nb]=i; fichier<<i<<", " ; cout<<i<<", " ;
652  nb++; if(nb==999){ cout<<"TO MUCH BAD CELLS"<<endl ; break;}}}
653  fichier<<"("<<nb<<")"<<endl; cout<<"("<<nb<<")"<<endl;}
654  fichier.close();
655 
656  if(file!="none"){
657  TFile::Open(file);
658  Int_t w=0 ;
659  Int_t c;
660  for(w=0; (w*9)<=nb; w++) {
661  if(9<=(nb-w*9)) c = 9 ;
662  else c = nb-9*w ;
663  Draw(bad, w*9, c) ;
664  }}
665 
666  nb=0;
667  }
668 
669 
670  //ANALYSIS
671  if (criterum < 3) TestCellEandN(pflag, Emin, Nsigma);
672  else if (criterum < 6)
673  TestCellShapes(pflag, Emin, Emax, Nsigma);
674 
675 
676  //RESULTS
677  if(criterum < 6) { nb=0;
678  cout<<"bad by lower value : "<<endl;
679  for(i=0;i<11520;i++) {
680  if(flag[i][criterum]==0 && exclu[i]==0){nb++;
681  cout<<i<<", " ;}} cout<<"("<<nb<<")"<<endl; nb=0;
682 
683  cout<<"bad by higher value : "<<endl;
684  for(i=0;i<11520;i++) {
685  if(flag[i][criterum]==2 && exclu[i]==0) {nb++;
686  cout<<i<<", " ;}} cout<<"("<<nb<<")"<<endl; nb=0;
687 
688  cout<<"total bad "<<endl;
689  for(i=0;i<11520;i++) {
690  if(flag[i][criterum]!=1 && exclu[i]==0) {
691  newBC[nb]=i;
692  nb++;
693  cout<<i<<", " ; }} cout<<"("<<nb<<")"<<endl;
694 
695 
696  //create a filtered file
697  KillCells(newBC,nb) ; nb=0;
698 
699  //write in a file the results
700  ofstream fichier(output, ios::out | ios::trunc);
701  if(fichier)
702  {
703  fichier <<"criterum : "<<criterum<<", Emin = "<<Emin<<" GeV"<<", Emax = "<<Emax<<" GeV"<<endl;
704  fichier<<"bad by lower value : "<<endl;
705  for(i=0;i<11520;i++) {
706  if(flag[i][criterum]==0 && exclu[i]==0){nb++;
707  fichier<<i<<", " ;}} fichier<<"("<<nb<<")"<<endl; nb=0;
708 
709  fichier<<"bad by higher value : "<<endl;
710  for(i=0;i<11520;i++) {
711  if(flag[i][criterum]==2 && exclu[i]==0) {nb++;
712  fichier<<i<<", " ;}} fichier<<"("<<nb<<")"<<endl; nb=0;
713 
714  fichier<<"total bad "<<endl;
715  for(i=0;i<11520;i++) {
716  if(flag[i][criterum]!=1 && exclu[i]==0) {
717  newBC[nb]=i;
718  nb++;
719  fichier<<i<<", " ; }} fichier<<"("<<nb<<")"<<endl;
720  fichier.close();
721  }
722  else
723  cerr << "opening error" << endl;
724  }
725 }
726 
727 
728 //_________________________________________________________________________
729 //_________________________________________________________________________
730 
737 void BCAnalysis(TString file, TString trigger = "default")
738 {
739  if(trigger=="default")
740  {
741  TFile::Open(file);
742  PeriodAnalysis(3, 8., 0.1, 2.5);
743  TFile::Open("filter.root");
744  PeriodAnalysis(2, 4., 0.1, 2.5);
745  TFile::Open("filter.root");
746  PeriodAnalysis(2, 4., 0.5, 2.5);
747  TFile::Open("filter.root");
748  PeriodAnalysis(1, 4., 0.1, 2.5);
749  TFile::Open("filter.root");
750  PeriodAnalysis(4, 4., 0.1,2.5);
751  }
752 
753  else
754  {
755  // you have the possibility to change analysis configuration in function of trigger type
756  TFile::Open(file);
757  PeriodAnalysis(3, 8., 0.1, 2.5);
758  TFile::Open("filter.root");
759  PeriodAnalysis(2, 4., 0.1, 2.5);
760  TFile::Open("filter.root");
761  PeriodAnalysis(2, 4., 0.5, 2.5);
762  TFile::Open("filter.root");
763  PeriodAnalysis(1, 4., 0.1, 2.5);
764  TFile::Open("filter.root");
765  PeriodAnalysis(4, 4., 0.1, 2.5);
766  }
767 
768  TFile::Open(file);
769  PeriodAnalysis(7,0.,0.,0.,file); //provide dead cells list from original file and draw bad cells candidate from indicated file
770 }
771 
772 //_________________________________________________________________________
773 //_________________________________________________________________________
774 
778 void BadChannelAnalysis(TString fCalorimeter = "EMCAL", TString period = "LHC11h",
779  TString pass = "pass1_HLT", TString trigger= "default")
780 {
781  Convert(fCalorimeter, period, pass, trigger);
782  TString inputfile = "/scratch/alicehp2/germain/QA/" + period;
783  inputfile += trigger;
784  inputfile += ".root";
785  BCAnalysis(inputfile, trigger);
786 }
void TestCellEandN(Int_t *pflag[11520][7], Double_t Emin=0.1, Double_t Nsigma=4., char *hname="hCellAmplitude", Int_t dnbins=200)
const char * title
Definition: MakeQAPdf.C:26
void Process(Int_t *pflag[11520][7], TH1 *inhisto, Double_t Nsigma=4., Int_t dnbins=200, Double_t dmaxval=-1.)
void TestCellShapes(Int_t *pflag[11520][7], Double_t fitEmin, Double_t fitEmax, Double_t Nsigma=4., char *hname="hCellAmplitude", Int_t dnbins=1000)
void PeriodAnalysis(Int_t criterum=7, Double_t Nsigma=4.0, Double_t Emin=0.1, Double_t Emax=1.0, TString file="none")
void BCAnalysis(TString file, TString trigger="default")
void BadChannelAnalysis(TString fCalorimeter="EMCAL", TString period="LHC11h", TString pass="pass1_HLT", TString trigger="default")
void Draw(Int_t cell[], Int_t iBC, Int_t nBC, const Int_t cellref=151)
void KillCells(Int_t filter[], Int_t nbc)
void Convert(TString fCalorimeter="EMCAL", TString period="LHC11h", TString pass="pass1_HLT", TString trigger="default")
TFile * file
void ExcludeCells(Int_t *pexclu[11520])