]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/TOF/trending/DrawTrendingTOFQA.C
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGPP / TOF / trending / DrawTrendingTOFQA.C
1 /*
2   fbellini@cern.ch - last update on 27/02/2014
3   Macro to draw the TOF QA trending plots by accessing the std tree.
4   To be mainly used with the automatic scripts to fill the QA repository.
5   Launch with 
6   aliroot -x -b -q "DrawTrendingTOFQA.C" 
7   The macro produces one png file for each trending variables
8   and a .root file with the histograms
9 */
10 Int_t DrawTrendingTOFQA(TString mergedTrendFile = "trending.root", // trending tree file 
11                         Bool_t displayAll = kFALSE) //set to kTRUE to display trending for expert plots
12 {
13   //
14   //reads merged trending.root file and draws trending plots from tree
15   //
16   if (!mergedTrendFile) {
17     Printf("Cannot open merged trend file with TOF QA");
18     return 1;
19   }
20   
21   char  outfilename[200]= "ProductionQA.hist.root";
22   TString plotDir(".");
23   // TString plotDir(Form("PlotsTrending"));
24   // gSystem->Exec(Form("mkdir %s",plotDir.Data()));
25
26   Int_t runNumber=0;
27   Double_t avTime=0., peakTime=0., spreadTime=0., peakTimeErr=0., spreadTimeErr=0.,negTimeRatio=0.,
28     avRawTime=0., peakRawTime=0., spreadRawTime=0., peakRawTimeErr=0., spreadRawTimeErr=0., 
29     avTot=0., peakTot=0.,spreadTot=0.,  peakTotErr=0.,spreadTotErr=0.,
30     orphansRatio=0., avL=0., negLratio=0.,
31     effPt1=0., effPt2=0., matchEffLinFit1Gev=0.,matchEffLinFit1GevErr=0.;
32   Double_t avDiffTime=0.,peakDiffTime=0., spreadDiffTime=0.,peakDiffTimeErr=0., spreadDiffTimeErr=0.,avT0fillRes=0.;
33    
34   Double_t avT0A=0.,peakT0A=0., spreadT0A=0.,peakT0AErr=0., spreadT0AErr=0.;
35   Double_t avT0C=0.,peakT0C=0., spreadT0C=0.,peakT0CErr=0., spreadT0CErr=0.;
36   Double_t avT0AC=0.,peakT0AC=0., spreadT0AC=0.,peakT0ACErr=0., spreadT0ACErr=0.;
37   Double_t avT0res=0.,peakT0res=0., spreadT0res=0.,peakT0resErr=0., spreadT0resErr=0.;
38   Float_t avMulti=0;
39   Float_t fractionEventsWHits=0.0;
40   Double_t goodChannelRatio=0.0;
41   
42   TFile * fin = TFile::Open(mergedTrendFile.Data());
43   TTree * ttree = (TTree*) fin->Get("trending");
44   if (!ttree){
45     Printf("Invalid trending tree.");
46     return 2;
47   }
48   ttree->SetBranchAddress("run",&runNumber);
49   ttree->SetBranchAddress("avMulti",&avMulti);
50   ttree->SetBranchAddress("goodChannelsRatio",&goodChannelRatio);   
51   ttree->SetBranchAddress("avTime",&avTime); //mean time
52   ttree->SetBranchAddress("peakTime",&peakTime); //main peak time after fit
53   ttree->SetBranchAddress("spreadTime",&spreadTime); //spread of main peak of time after fit
54   ttree->SetBranchAddress("peakTimeErr",&peakTimeErr); //main peak time after fit error
55   ttree->SetBranchAddress("spreadTimeErr",&spreadTimeErr); //spread of main peak of time after fit error
56   ttree->SetBranchAddress("negTimeRatio",&negTimeRatio); //negative time ratio
57   ttree->SetBranchAddress("avRawTime",&avRawTime); //mean raw time
58   ttree->SetBranchAddress("peakRawTime",&peakRawTime); //mean peak of raw time after fit
59   ttree->SetBranchAddress("spreadRawTime",&spreadRawTime); //spread of main peak of raw time after fit
60   ttree->SetBranchAddress("peakRawTimeErr",&peakRawTimeErr); //main peak raw  time after fit error
61   ttree->SetBranchAddress("spreadRawTimeErr",&spreadRawTimeErr); //spread of  raw main peak of time after fit error
62   ttree->SetBranchAddress("avTot",&avTot); //main peak tot
63   ttree->SetBranchAddress("peakTot",&peakTot); // main peak of tot after fit
64   ttree->SetBranchAddress("spreadTot",&spreadTot); //spread of main peak of tot after fit
65   ttree->SetBranchAddress("peakTotErr",&peakTotErr); // main peak of tot after fit
66   ttree->SetBranchAddress("spreadTotErr",&spreadTotErr); //spread of main peak of tot after fit
67   ttree->SetBranchAddress("orphansRatio",&orphansRatio); //orphans ratio
68   ttree->SetBranchAddress("avL",&avL); //mean track length
69   ttree->SetBranchAddress("negLratio",&negLratio);//ratio of tracks with track length <350 cm
70   ttree->SetBranchAddress("effPt1",&effPt1);//matching eff at 1 GeV/c
71   ttree->SetBranchAddress("effPt2",&effPt2); //matching eff at 2 GeV/c
72   ttree->SetBranchAddress("matchEffLinFit1Gev",&matchEffLinFit1Gev);//matching eff fit param 
73   ttree->SetBranchAddress("matchEffLinFit1GevErr",&matchEffLinFit1GevErr);////matching eff fit param error
74   ttree->SetBranchAddress("avPiDiffTime",&avDiffTime); //mean t-texp
75   ttree->SetBranchAddress("peakPiDiffTime",&peakDiffTime); //main peak t-texp after fit
76   ttree->SetBranchAddress("spreadPiDiffTime",&spreadDiffTime); //spread of main peak t-texp after fit
77   ttree->SetBranchAddress("peakPiDiffTimeErr",&peakDiffTimeErr); //main peak t-texp after fit error
78   ttree->SetBranchAddress("spreadPiDiffTimeErr",&spreadDiffTimeErr); //spread of main peak of t-texp after fit error
79   ttree->SetBranchAddress("avT0A",&avT0A); //main peak t0A
80   ttree->SetBranchAddress("peakT0A",&peakT0A); // main peak of t0A after fit
81   ttree->SetBranchAddress("spreadT0A",&spreadT0A); //spread of main peak of t0A after fit
82   ttree->SetBranchAddress("peakT0AErr",&peakT0AErr); // main peak of t0A after fit
83   ttree->SetBranchAddress("spreadT0AErr",&spreadT0AErr); //spread of main peak of t0A after fit
84   ttree->SetBranchAddress("avT0C",&avT0C); //main peak t0C
85   ttree->SetBranchAddress("peakT0C",&peakT0C); // main peak of t0C after fit
86   ttree->SetBranchAddress("spreadT0C",&spreadT0C); //spread of main peak of t0C after fit
87   ttree->SetBranchAddress("peakT0CErr",&peakT0CErr); // main peak of t0C after fit
88   ttree->SetBranchAddress("spreadT0CErr",&spreadT0CErr); //spread of main peak of t0C after fit
89   ttree->SetBranchAddress("avT0AC",&avT0AC); //main peak t0AC
90   ttree->SetBranchAddress("peakT0AC",&peakT0AC); // main peak of t0AC after fit
91   ttree->SetBranchAddress("spreadT0AC",&spreadT0AC); //spread of main peak of t0AC after fit
92   ttree->SetBranchAddress("peakT0ACErr",&peakT0ACErr); // main peak of t0AC after fit
93   ttree->SetBranchAddress("spreadT0ACErr",&spreadT0ACErr); //spread of main peak of t0AC after fit
94   ttree->SetBranchAddress("avT0res",&avT0res); //main peak t0AC
95   ttree->SetBranchAddress("peakT0res",&peakT0res); // main peak of t0AC after fit
96   ttree->SetBranchAddress("spreadT0res",&spreadT0res); //spread of main peak of t0AC after fit
97   ttree->SetBranchAddress("peakT0resErr",&peakT0resErr); // main peak of t0AC after fit
98   ttree->SetBranchAddress("spreadT0resErr",&spreadT0resErr); //spread of main peak of t0AC after fit
99   ttree->SetBranchAddress("avT0fillRes",&avT0fillRes); //t0 fill res
100
101   Int_t nRuns=ttree->GetEntries();
102   TList lista;
103    
104   TH1F * hAvMulti=new TH1F("hAvMulti","Average multiplicity of matched tracks <N_{TOF}>;; <N_{TOF}>", nRuns,0., nRuns);//, 600, 0. , 600.);
105   hAvMulti->SetDrawOption("E1");
106   hAvMulti->SetMarkerStyle(20);
107   hAvMulti->SetMarkerColor(kBlue);
108         
109   TH1F * hAvDiffTimeVsRun=new TH1F("hAvDiffTimeVsRun","Mean t-t_{exp} (no fit);run;<t^{TOF}-t_{exp,#pi}> (ps)",nRuns,0., nRuns);//, 600, 0. , 600.);
110   hAvDiffTimeVsRun->SetDrawOption("E1");
111   hAvDiffTimeVsRun->SetMarkerStyle(20);
112   hAvDiffTimeVsRun->SetMarkerColor(kBlue);
113   //   hAvTimeVsRun->GetYaxis()->SetRangeUser(0.0, 50.0);
114
115   TH1F * hPeakDiffTimeVsRun=new TH1F("hPeakDiffTimeVsRun","t-t_{exp} (gaussian fit) ;; <t^{TOF}-t_{exp,#pi}> (ps)",nRuns,0., nRuns);//,600, 0. , 600. );
116   hPeakDiffTimeVsRun->SetDrawOption("E1");
117   hPeakDiffTimeVsRun->SetMarkerStyle(20);
118   hPeakDiffTimeVsRun->SetMarkerColor(kBlue);
119    
120   TH1F * hSpreadDiffTimeVsRun=new TH1F("hSpreadDiffTimeVsRun","#sigma(t-t_{exp}) (gaussian fit);; #sigma(t^{TOF}-t_{exp,#pi}) (ns)",nRuns,0., nRuns);//, 100, 0. , 100.);
121   hSpreadDiffTimeVsRun->SetDrawOption("E1");
122   hSpreadDiffTimeVsRun->SetMarkerStyle(20);
123   hSpreadDiffTimeVsRun->SetMarkerColor(kBlue);
124
125   TH1F * hAvTimeVsRun=new TH1F("hAvTimeVsRun","<t^{TOF}>;;<t^{TOF}> (ns)",nRuns,0., nRuns);//, 600, 0. , 600.);
126   hAvTimeVsRun->SetDrawOption("E1");
127   hAvTimeVsRun->SetMarkerStyle(20);
128   hAvTimeVsRun->SetMarkerColor(kBlue);
129   //   hAvTimeVsRun->GetYaxis()->SetRangeUser(0.0, 50.0);
130
131   TH1F * hPeakTimeVsRun=new TH1F("hPeakTimeVsRun","Peak value of t^{TOF} (landau fit);;t_{peak}^{TOF} (ns)",nRuns,0., nRuns);//,600, 0. , 600. );
132   hPeakTimeVsRun->SetDrawOption("E1");
133   hPeakTimeVsRun->SetMarkerStyle(20);
134   hPeakTimeVsRun->SetMarkerColor(kBlue);
135    
136   TH1F * hSpreadTimeVsRun=new TH1F("hSpreadTimeVsRun","Spread of t^{TOF} (landau fit);; #sigma(t^{TOF}) (ns)",nRuns,0., nRuns);//, 100, 0. , 100.);
137   hSpreadTimeVsRun->SetDrawOption("E1");
138   hSpreadTimeVsRun->SetMarkerStyle(20);
139   hSpreadTimeVsRun->SetMarkerColor(kBlue);
140   
141   TH1F * hAvRawTimeVsRun=new TH1F("hAvRawTimeVsRun","Peak value of raw t^{TOF};;<t_{raw}^{TOF}> (ns)",nRuns,0., nRuns);//, 600, 0. , 600.);
142   hAvRawTimeVsRun->SetDrawOption("E1");
143   hAvRawTimeVsRun->SetMarkerStyle(21);
144   hAvRawTimeVsRun->SetMarkerColor(kGreen);
145
146   TH1F * hPeakRawTimeVsRun=new TH1F("hPeakRawTimeVsRun","Peak value of raw t^{TOF} (landau fit);;t_{peak,raw}^{TOF} (ns)",nRuns,0., nRuns);//, 600, 0. , 600.);
147   hPeakRawTimeVsRun->SetDrawOption("E1");
148   hPeakRawTimeVsRun->SetMarkerStyle(21);
149   hPeakRawTimeVsRun->SetMarkerColor(kGreen);
150
151   TH1F * hSpreadRawTimeVsRun=new TH1F("hSpreadRawTimeVsRun","Spread of raw t^{TOF} (landau fit);;#sigma(t_{raw}^{TOF}) (ns)",nRuns,0., nRuns);//, 100, 0. , 100.);
152   hSpreadRawTimeVsRun->SetDrawOption("E1");
153   hSpreadRawTimeVsRun->SetMarkerStyle(21);
154   hSpreadRawTimeVsRun->SetMarkerColor(kGreen);
155    
156   TH1F * hAvTotVsRun=new TH1F("hAvTotVsRun","<ToT> (no fit);run;<ToT> (ns)",nRuns,0., nRuns);//, 50, 0. , 50.);
157   hAvTotVsRun->SetDrawOption("E1");
158   hAvTotVsRun->SetMarkerStyle(22);
159    
160   TH1F * hPeakTotVsRun=new TH1F("hPeakTotVsRun","<ToT> (gaussian fit);;ToT_{peak} (ns)",nRuns,0., nRuns);//, 50, 0. , 50.);
161   hPeakTotVsRun->SetDrawOption("E1");
162   hPeakTotVsRun->SetMarkerStyle(22);
163    
164   TH1F * hSpreadTotVsRun=new TH1F("hSpreadTotVsRun","#sigma(ToT) (gaussian fit);#sigma(ToT) (ns)",nRuns,0., nRuns);//, 50, 0. , 50.);
165   hSpreadTotVsRun->SetDrawOption("E1");
166   hSpreadTotVsRun->SetMarkerStyle(22);
167    
168   TH1F * hNegTimeRatioVsRun=new TH1F("hNegTimeRatioVsRun","Ratio of tracks with t^{TOF}<12.5 ns; ; ratio of tracks with t^{TOF}<12.5 ns (%)",nRuns, 0., nRuns);//, 100, 0. , 100.);
169   hNegTimeRatioVsRun->SetDrawOption("E");
170
171   TH1F * hOrphansRatioVsRun=new TH1F("hOrphansRatioVsRun","Ratio of orphans (hits with ToT=0); ; ratio of orphans (%)",nRuns, 0., nRuns);//, 1000, 0. , 100.);
172   hOrphansRatioVsRun->SetDrawOption("E");
173
174   TH1F * hMeanLVsRun=new TH1F("hMeanLVsRun","Average track length;; <L> (cm)",nRuns, 0., nRuns);//, 350, 350. , 700.);
175   hMeanLVsRun->SetDrawOption("E");
176   TH1F * hNegLRatioVsRun=new TH1F("hNegLRatioVsRun","Ratio of tracks with L<350 cm;; ratio of tracks with L<350 cm (%)",nRuns, 0., nRuns);//, 1000, 0. , 100.);
177   hNegLRatioVsRun->SetDrawOption("E");
178   TH1F * hMatchEffVsRun=new TH1F("hMatchEffVsRun","Matching efficiency (linear fit for p_{T}>1.0 GeV/c);;matching efficiency (pT>1.0 GeV/c)",nRuns, 0., nRuns);//, 100, 0. , 1.);
179   hMatchEffVsRun->SetDrawOption("E");
180   TH1F * hMatchEffVsRunNormToGoodCh=new TH1F("hMatchEffVsRunNormToGoodCh","Matching efficiency normalized to percentage of TOF good channels;;matching efficiency (pT>1.0 GeV/c)",nRuns, 0., nRuns);//, 100, 0. , 1.);
181   hMatchEffVsRunNormToGoodCh->SetDrawOption("E");
182         
183   TH1F * hMatchEffVsRun1=new TH1F("hMatchEffVsRun1","Matching efficiency (value for p_{T}=1.0 GeV/c);;matching efficiency (pT=1.0 GeV/c)",nRuns, 0., nRuns);
184   hMatchEffVsRun1->SetDrawOption("E");
185   TH1F * hPeakT0AVsRun=new TH1F("hPeakT0AVsRun","Peak value of T0A (gaussian fit);;t0A (ps)",nRuns,0., nRuns);
186   TH1F * hPeakT0CVsRun=new TH1F("hPeakT0CVsRun","Peak value of T0C (gaussian fit);;t0AC (ps)",nRuns,0., nRuns);
187   TH1F * hPeakT0ACVsRun=new TH1F("hPeakT0ACVsRun","Peak value of T0AC (gaussian fit);;t0AC (ps)",nRuns,0., nRuns);
188   TH1F * hT0fillResVsRun=new TH1F("hT0fillResVsRun","t0_fill spread;;t0_spread (ps)",nRuns,0., nRuns);
189         
190   TH1F * hGoodChannelsRatio=new TH1F("hGoodChannelsRatio","Fraction of TOF good channels;;fraction of good channels",nRuns, 0., nRuns);//, 100, 0. , 1.);
191   hGoodChannelsRatio->SetDrawOption("E");
192         
193   lista.Add(hAvMulti);
194   lista.Add(hAvDiffTimeVsRun);
195   lista.Add(hPeakDiffTimeVsRun);
196   lista.Add(hSpreadDiffTimeVsRun);
197   lista.Add(hAvTimeVsRun);
198   lista.Add(hPeakTimeVsRun);
199   lista.Add(hSpreadTimeVsRun);
200   lista.Add(  hAvRawTimeVsRun);
201   lista.Add(  hPeakRawTimeVsRun);
202   lista.Add(  hSpreadRawTimeVsRun); 
203   lista.Add(  hAvTotVsRun);
204   lista.Add(  hPeakTotVsRun);
205   lista.Add(  hSpreadTotVsRun);
206   lista.Add(  hNegTimeRatioVsRun);
207   lista.Add(  hOrphansRatioVsRun);
208   lista.Add( hMeanLVsRun);
209   lista.Add(  hNegLRatioVsRun);
210   lista.Add(  hMatchEffVsRun);
211   lista.Add(hMatchEffVsRunNormToGoodCh);
212   lista.Add(hPeakT0AVsRun);
213   lista.Add(hPeakT0CVsRun);
214   lista.Add(hPeakT0ACVsRun);
215   lista.Add(hT0fillResVsRun);
216   lista.Add(hGoodChannelsRatio);
217
218   char runlabel[6];
219    
220   for (Int_t irun=0;irun<nRuns;irun++){
221     ttree->GetEntry(irun);
222     
223     sprintf(runlabel,"%i",runNumber);
224     
225     hAvMulti->SetBinContent(irun+1, avMulti);
226     hAvMulti->GetXaxis()->SetBinLabel(irun+1,runlabel);
227
228     hAvDiffTimeVsRun->SetBinContent(irun+1, avDiffTime);
229     hAvDiffTimeVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
230
231     hPeakDiffTimeVsRun->SetBinContent(irun+1,peakDiffTime);
232     hPeakDiffTimeVsRun->SetBinError(irun+1,peakDiffTimeErr);
233     hPeakDiffTimeVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
234
235     hSpreadDiffTimeVsRun->SetBinContent(irun+1,spreadDiffTime);
236     hSpreadDiffTimeVsRun->SetBinError(irun+1,spreadDiffTimeErr);
237     hSpreadDiffTimeVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
238
239     hAvTimeVsRun->SetBinContent(irun+1, avTime);
240     hAvTimeVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
241
242     hPeakTimeVsRun->SetBinContent(irun+1,peakTime);
243     hPeakTimeVsRun->SetBinError(irun+1,peakTimeErr);
244     hPeakTimeVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
245
246     hSpreadTimeVsRun->SetBinContent(irun+1,spreadTime);
247     hSpreadTimeVsRun->SetBinError(irun+1,spreadTimeErr);
248     hSpreadTimeVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
249
250     hAvRawTimeVsRun->SetBinContent(irun+1, avRawTime);
251     hAvRawTimeVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
252     
253     hPeakRawTimeVsRun->SetBinContent(irun+1,peakRawTime);
254     hPeakRawTimeVsRun->SetBinError(irun+1,peakRawTimeErr);
255     hPeakRawTimeVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
256
257     hSpreadRawTimeVsRun->SetBinContent(irun+1,spreadRawTime);
258     hSpreadRawTimeVsRun->SetBinError(irun+1,spreadRawTimeErr);
259     hSpreadRawTimeVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
260
261     hAvTotVsRun->SetBinContent(irun+1,avTot);
262     hAvTotVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
263
264     hPeakTotVsRun->SetBinContent(irun+1,peakTot);
265     hPeakTotVsRun->SetBinError(irun+1,peakTotErr);
266     hPeakTotVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
267
268     hSpreadTotVsRun->SetBinContent(irun+1,spreadTot);
269     hSpreadTotVsRun->SetBinError(irun+1,spreadTotErr);
270     hSpreadTotVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
271     
272     hNegTimeRatioVsRun->SetBinContent(irun+1,negTimeRatio);
273     hNegTimeRatioVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
274     
275     hOrphansRatioVsRun->SetBinContent(irun+1,orphansRatio);
276     hOrphansRatioVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
277     
278     hMeanLVsRun->SetBinContent(irun+1,avL);
279     hMeanLVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
280     
281     hNegLRatioVsRun->SetBinContent(irun+1,negLratio);
282     hNegLRatioVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
283
284     hMatchEffVsRun->SetBinContent(irun+1,matchEffLinFit1Gev);
285     hMatchEffVsRun->SetBinError(irun+1,matchEffLinFit1GevErr);
286     hMatchEffVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
287     hMatchEffVsRun->SetLineColor(kRed);
288     hMatchEffVsRun->SetLineWidth(2);
289
290     if (goodChannelRatio>0)
291       hMatchEffVsRunNormToGoodCh->SetBinContent(irun+1,matchEffLinFit1Gev/goodChannelRatio);
292     else 
293       hMatchEffVsRunNormToGoodCh->SetBinContent(irun+1, 0.0);
294     hMatchEffVsRunNormToGoodCh->SetBinError(irun+1,matchEffLinFit1GevErr);
295     hMatchEffVsRunNormToGoodCh->GetXaxis()->SetBinLabel(irun+1,runlabel);
296     hMatchEffVsRunNormToGoodCh->SetLineColor(kBlue);
297     hMatchEffVsRunNormToGoodCh->SetLineWidth(2);
298      
299     hGoodChannelsRatio->SetBinContent(irun+1, goodChannelRatio);
300     hGoodChannelsRatio->SetLineColor(kMagenta+2);
301     hGoodChannelsRatio->SetLineWidth(2);
302     hGoodChannelsRatio->GetXaxis()->SetBinLabel(irun+1,runlabel);
303
304     hPeakT0AVsRun->SetBinContent(irun+1,peakT0A);
305     hPeakT0AVsRun->SetBinError(irun+1,spreadT0A);
306     hPeakT0AVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
307
308     hPeakT0CVsRun->SetBinContent(irun+1,peakT0C);
309     hPeakT0CVsRun->SetBinError(irun+1,spreadT0C);
310     hPeakT0CVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
311     
312     hPeakT0ACVsRun->SetBinContent(irun+1,peakT0AC);
313     hPeakT0ACVsRun->SetBinError(irun+1,spreadT0AC);
314     hPeakT0ACVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
315
316     hT0fillResVsRun->SetBinContent(irun+1,avT0fillRes);
317     hT0fillResVsRun->GetXaxis()->SetBinLabel(irun+1,runlabel);
318   }
319   
320   TFile * fout=new TFile(outfilename,"recreate");
321   fout->cd();
322   lista.Write();
323   fout->Close();
324     
325   gStyle->SetOptStat(10);
326
327   TCanvas* cPeakDiffTimeVsRun = new TCanvas("cPeakDiffTimeVsRun","cPeakDiffTimeVsRun", 50,50,1050, 550);
328   hPeakDiffTimeVsRun->GetYaxis()->SetRangeUser(-50.,50.);
329   hPeakDiffTimeVsRun->Draw();
330   cPeakDiffTimeVsRun->Print(Form("%s/cPeakDiffTimeVsRun.png",plotDir.Data()));
331         
332   TCanvas* cSpreadDiffTimeVsRun = new TCanvas("cSpreadDiffTimeVsRun","cSpreadDiffTimeVsRun", 50,50,1050, 550);
333   hSpreadDiffTimeVsRun->GetYaxis()->SetRangeUser(200.,400.);
334   hSpreadDiffTimeVsRun->Draw();
335   cSpreadDiffTimeVsRun->Print(Form("%s/cSpreadDiffTimeVsRun.png",plotDir.Data()));
336   
337   TCanvas* cMatchEffVsRun = new TCanvas("cMatchEffVsRun","cMatchEffVsRun",50, 50, 1050, 550);
338   hMatchEffVsRun->GetYaxis()->SetRangeUser(0.,1.);
339   hMatchEffVsRun->Draw();
340   cMatchEffVsRun->Print(Form("%s/cMatchEffVsRun.png",plotDir.Data()));
341   
342   TCanvas* cMatchEffNormToGoodCh = new TCanvas("cMatchEffNormToGoodCh","cMatchEffNormToGoodCh",50, 50,1050, 550);
343   hMatchEffVsRunNormToGoodCh->GetYaxis()->SetRangeUser(0.,1.);
344   hMatchEffVsRunNormToGoodCh->Draw();
345   cMatchEffNormToGoodCh->Print(Form("%s/cMatchEffNormToGoodCh.png",plotDir.Data()));
346
347   TCanvas* cGoodCh = new TCanvas("cGoodCh","cGoodCh",50, 50,1050, 550);
348   hGoodChannelsRatio->GetYaxis()->SetRangeUser(0.75,1.);
349   hGoodChannelsRatio->Draw();
350   cGoodCh->Print(Form("%s/cGoodCh.png",plotDir.Data()));
351
352   if (displayAll) {     
353     TCanvas* cPeakT0AVsRun = new TCanvas("cPeakT0AVsRun","cPeakT0AVsRun", 50,50,1050, 550);
354     hPeakT0AVsRun->Draw();
355     cPeakT0AVsRun->Print(Form("%s/cPeakT0AVsRun.png",plotDir.Data()));
356   
357     TCanvas* cPeakT0CVsRun = new TCanvas("cPeakT0CVsRun","cPeakT0CVsRun", 50,50,1050, 550);
358     hPeakT0CVsRun->Draw();
359     cPeakT0CVsRun->Print(Form("%s/cPeakT0CVsRun.png",plotDir.Data()));
360   
361     TCanvas* cPeakT0ACVsRun = new TCanvas("cPeakT0ACVsRun","cPeakT0ACVsRun", 50,50,1050, 550);
362     hPeakT0ACVsRun->Draw();
363     cPeakT0ACVsRun->Print(Form("%s/cPeakT0ACVsRun.png",plotDir.Data()));
364   
365     TCanvas* cT0fillResVsRun = new TCanvas("cT0fillResVsRun","cT0fillResVsRun", 50,50,1050, 550);
366     hT0fillResVsRun->Draw();
367     cT0fillResVsRun->Print(Form("%s/cT0fillResVsRun.png",plotDir.Data()));
368
369     TCanvas* cAvDiffTimeVsRun = new TCanvas("cAvDiffTimeVsRun","cAvDiffTimeVsRun",50,50,1050, 550);
370     gPad->SetGridx();
371     gPad->SetGridy();
372     hAvDiffTimeVsRun->Draw();
373     cAvDiffTimeVsRun->Print(Form("%s/cAvDiffTimeVsRun.png",plotDir.Data()));
374           
375     TCanvas* cAvTimeVsRun = new TCanvas("cAvTimeVsRun","cAvTimeVsRun", 50,50,1050, 550);
376     hAvTimeVsRun->Draw();
377     cAvTimeVsRun->Print(Form("%s/cAvTimeVsRun.png",plotDir.Data()));
378           
379     TCanvas* cPeakTimeVsRun = new TCanvas("cPeakTimeVsRun","cPeakTimeVsRun", 50,50,1050, 550);
380     hPeakTimeVsRun->Draw();
381     cPeakTimeVsRun->Print(Form("%s/cPeakTimeVsRun.png",plotDir.Data()));
382           
383     TCanvas* cSpreadTimeVsRun = new TCanvas("cSpreadTimeVsRun","cSpreadTimeVsRun", 50,50,1050, 550);
384     hSpreadTimeVsRun->Draw();
385     cSpreadTimeVsRun->Print(Form("%s/cSpreadTimeVsRun.png",plotDir.Data()));
386           
387     TCanvas* cAvRawTimeVsRun = new TCanvas("cAvRawTimeVsRun","cAvRawTimeVsRun", 50,50,1050, 550);
388     hAvRawTimeVsRun->Draw();
389     cAvRawTimeVsRun->Print(Form("%s/cAvRawTimeVsRun.png",plotDir.Data()));
390           
391     TCanvas* cPeakRawTimeVsRun = new TCanvas("cPeakRawTimeVsRun","cPeakRawTimeVsRun", 50,50,1050, 550);
392     hPeakRawTimeVsRun->Draw();
393     cPeakRawTimeVsRun->Print(Form("%s/cPeakRawTimeVsRun.png",plotDir.Data()));
394           
395     TCanvas* cSpreadRawTimeVsRun = new TCanvas("cSpreadRawTimeVsRun","cSpreadRawTimeVsRun", 50,50,1050, 550);
396     hSpreadRawTimeVsRun->Draw();
397     cSpreadRawTimeVsRun->Print(Form("%s/cSpreadRawTimeVsRun.png",plotDir.Data()));
398           
399     TCanvas* cAvTotVsRun = new TCanvas("cAvTotVsRun","cAvTotVsRun", 50,50,1050, 550);
400     hAvTotVsRun->Draw();
401     cAvTotVsRun->Print(Form("%s/cAvTotVsRun.png",plotDir.Data()));
402           
403     TCanvas* cPeakTotVsRun = new TCanvas("cPeakTotVsRun","cPeakTotVsRun", 50,50,1050, 550);
404     hPeakTotVsRun->Draw();
405     cPeakTotVsRun->Print(Form("%s/cPeakTotVsRun.png",plotDir.Data()));
406           
407     TCanvas* cSpreadTotVsRun = new TCanvas("cSpreadTotVsRun","cSpreadTotVsRun", 50,50,1050, 550);
408     hSpreadTotVsRun->Draw();
409     cSpreadTotVsRun->Print(Form("%s/cSpreadTotVsRun.png",plotDir.Data()));
410           
411     TCanvas* cNegTimeRatioVsRun = new TCanvas("cNegTimeRatioVsRun","cNegTimeRatioVsRun", 50,50,1050, 550);
412     hNegTimeRatioVsRun->Draw();
413     cNegTimeRatioVsRun->Print(Form("%s/cNegTimeRatioVsRun.png",plotDir.Data()));
414           
415     TCanvas* cOrphansRatioVsRun = new TCanvas("cOrphansRatioVsRun","cOrphansRatioVsRun", 50,50,1050, 550);
416     hOrphansRatioVsRun->Draw();
417     cOrphansRatioVsRun->Print(Form("%s/cOrphansRatioVsRun.png",plotDir.Data()));
418           
419     TCanvas* cMeanLVsRun = new TCanvas("cMeanLVsRun","cMeanLVsRun", 50,50,1050, 550);
420     hMeanLVsRun->Draw();
421     cMeanLVsRun->Print(Form("%s/cMeanLVsRun.png",plotDir.Data()));
422           
423     TCanvas* cNegLRatioVsRun = new TCanvas("cNegLRatioVsRun","cNegLRatioVsRun", 50,50,1050, 550);
424     hNegLRatioVsRun->Draw();
425     cNegLRatioVsRun->Print(Form("%s/cNegLRatioVsRun.png",plotDir.Data()));
426   }
427         
428   return 0;
429 }