]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ZDC/ZDCQAtrending.C
Updating Shuttle preprocessor
[u/mrichter/AliRoot.git] / ZDC / ZDCQAtrending.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2
3 #include <stdio.h>
4 #include <stdlib.h>
5 #include <TROOT.h>
6 #include <Riostream.h>
7 #include <TClassTable.h>
8 #include <TStyle.h>
9 #include <TMath.h>
10 #include <TFile.h>
11 #include <TCanvas.h>
12 #include <TH1.h>
13 #include <TH2.h>
14 #include <TProfile.h>
15 #include <TLine.h>
16 #include <TGrid.h>
17 #include <TBits.h>
18 #include <TChain.h>
19 #include <TNtuple.h>
20 #include <TTree.h>
21 #include <TBranch.h>
22 #include <TFileMerger.h>
23 #include <TGridResult.h>
24 #include <TSystem.h>
25 #include <TLegend.h>
26
27 #endif
28
29 void MakePlots(TString ntupleFileName);
30
31 void ZDCQAtrending(TString period,
32                    TString recoPass,
33                    TString qaTrain = "QA",
34                    Bool_t useOnlyMerged = kTRUE,
35                    Int_t firstRun = 0,
36                    Int_t lastRun = 999999999,
37                    TString runListFile = "",
38                    TString fileName = "QAresults.root"){
39
40   gStyle->SetOptStat(0);
41
42   TGrid::Connect("alien:");
43   Int_t year = 0;
44   if(period.Contains("LHC09")) year = 2009;
45   else if(period.Contains("LHC10")) year = 2010;
46   else if(period.Contains("LHC11")) year = 2011;
47   
48   Bool_t useExternalList = kFALSE;
49   Int_t runList[10000];
50   Int_t totRuns = 0;
51   if(runListFile.Length()>0){
52     if(!gSystem->Exec(Form("ls -l %s > /dev/null 2>&1",runListFile.Data()))){
53       printf("Use Run List from %s  --- runs to be analyzed:\n",runListFile.Data());
54       useExternalList = kTRUE;
55       FILE* rfil = fopen(runListFile.Data(),"r");
56       Int_t nrun;
57       while(!feof(rfil)){
58         int stat = fscanf(rfil,"%d, ",&nrun);
59         if(feof(rfil)) break;
60         runList[totRuns++] = nrun;
61       }
62       for(Int_t ir = 0; ir<totRuns; ir++){
63         printf("%d\n",runList[ir]);
64       }
65     }
66     else{
67       printf("File with run list does not exist\n");
68     }
69   }
70
71   TString outFilNam = Form("ZDCtrend_%s_%s_%s.root",period.Data(),recoPass.Data(),qaTrain.Data());
72
73  
74   const Int_t nVariables = 27;
75   TNtuple* ntzdc = new TNtuple("ntzdc","ZDC trending",
76   "nrun:meanZNC:meanZPC:meanZNA:meanZPA:meanZEM1:meanZEM2:emeanZNC:emeanZPC:emeanZNA:emeanZPA:emeanZEM1:emeanZEM2:"
77   "ZNClg:ZNAlg:eZNClg:eZNAlg:pmcZNClg:pmcZNAlg:epmcZNClg:epmcZNAlg:xZNC:yZNC:xZNA:yZNA:tdcSum:tdcDiff");
78   Float_t xnt[nVariables];
79
80   TBits* readRun = new TBits(999999);
81   readRun->ResetAllBits();
82   if(!useExternalList){
83     if(!gSystem->Exec(Form("ls -l %s > /dev/null 2>&1",outFilNam.Data()))){
84       TFile* oldfil = new TFile(outFilNam.Data());
85       TNtuple* ntmp = (TNtuple*)oldfil->Get("ntzdc");
86       Bool_t isOK = kFALSE;
87       if(ntmp){
88         if(ntmp->GetNvar() == ntzdc->GetNvar()){
89           isOK = kTRUE;
90           TObjArray* arr1 = (TObjArray*)ntzdc->GetListOfBranches();
91           TObjArray* arr2 = (TObjArray*)ntmp->GetListOfBranches();
92           for(Int_t iV = 0; iV<ntmp->GetNvar(); iV++){
93             TString vnam1 = arr1->At(iV)->GetName();
94             TString vnam2 = arr2->At(iV)->GetName();
95             if(vnam1 != vnam2) isOK = kFALSE;
96             ntmp->SetBranchAddress(vnam2.Data(),&xnt[iV]);
97           }
98           if(isOK){
99             for(Int_t nE = 0; nE<ntmp->GetEntries(); nE++){
100               ntmp->GetEvent(nE);
101               Int_t theRun = (Int_t)(xnt[0]+0.0001);
102               readRun->SetBitNumber(theRun);
103               ntzdc->Fill(xnt);
104             }
105           }
106         }
107       }
108       if(!isOK){
109         printf("Ntuple in local file not OK -> will be recreated\n");
110       }
111       oldfil->Close();
112       delete oldfil;
113     }
114   }
115
116   if(!gGrid||!gGrid->IsConnected()) {
117     printf("gGrid not found! exit macro\n");
118     return;
119   }
120
121   TString  path = Form("/alice/data/%d/%s/",year,period.Data());
122   TGridResult *gr  =  gGrid->Query(path,fileName);
123   Int_t nFiles  =  gr->GetEntries();
124   printf(" --->%d merged QAresults.root files found in %s\n", nFiles,path.Data());
125   if(nFiles < 1) return;
126
127   Int_t nAnalyzedFiles = 0;
128   if(nFiles > 1){
129     for(Int_t iFil  =  0; iFil <nFiles ; iFil++) { 
130       TString fileNameLong = Form("%s",gr->GetKey(iFil,"turl"));
131       if(!fileNameLong.Contains(recoPass.Data())) continue;
132 // Commented Sep. 2011
133 //      if(!fileNameLong.Contains(qaTrain.Data())) continue;
134 //      if(fileNameLong.Contains("TRD") || fileNameLong.Contains("EMCAL")) continue;
135       TString runNumber = fileNameLong;
136       runNumber.ReplaceAll(Form("alien:///alice/data/%d/%s/",year,period.Data()),"");
137       runNumber.Remove(9,runNumber.Sizeof());
138    
139       Int_t iRun = atoi(runNumber.Data());
140       //printf("  runNumber: %s -> iRun %d\n",runNumber.Data(), iRun);
141       if(useExternalList){
142         Bool_t keepRun = kFALSE;
143         for(Int_t ir = 0; ir<totRuns; ir++){
144           if(iRun == runList[ir]){
145             keepRun = kTRUE;
146             break;
147           }
148         }
149         if(!keepRun) continue;
150       }
151       if(readRun->TestBitNumber(iRun)){ 
152         printf("Run %d already in local ntuple -> skipping it\n",iRun);
153         continue;
154       }
155       if(iRun<firstRun) continue;
156       if(iRun>lastRun) continue;    
157
158 printf("  fileNameLog: %s\n",fileNameLong.Data());
159       if(useOnlyMerged){
160         TString isMerged = fileNameLong;
161         //
162         if(!isMerged.Contains(recoPass.Data())) continue;
163         if(isMerged.Contains("Stage")) continue;
164 printf("  ->  isMerged: %s\n",isMerged.Data());
165         isMerged.Remove(0,isMerged.Sizeof()-20); 
166 printf("  ->  isMerged: %s\n",isMerged.Data());
167         if(!isMerged.Contains(qaTrain) || !isMerged.Contains("QAresults.root")) continue;
168 printf("isMerged %s\n",isMerged.Data());
169         
170       }
171       printf("Open File %s  \n",fileNameLong.Data());
172       
173
174       TFile* f = TFile::Open(fileNameLong.Data());  
175
176       TDirectoryFile* df = (TDirectoryFile*)f->Get("ZDC_Performance");
177       if(!df){
178         printf("Run %d ZDC_Performance MISSING -> Exit\n",iRun);
179         continue;
180       }
181       TList* l = (TList*)df->Get("QAZDCHists");
182       if(!df){
183         printf("Run %d QAZDCHists TList MISSING -> Exit\n",iRun);
184         continue;
185       }  
186
187       nAnalyzedFiles++;
188       if(!useOnlyMerged) readRun->SetBitNumber(iRun);
189
190       TH1F *fhTDCZNSum     = (TH1F*)l->FindObject("fhTDCZNSum");     
191       TH1F *fhTDCZNDiff    = (TH1F*)l->FindObject("fhTDCZNDiff");    
192       TH1F *fhZNCSpectrum  = (TH1F*)l->FindObject("fhZNCSpectrum");  
193       TH1F *fhZNASpectrum  = (TH1F*)l->FindObject("fhZNASpectrum");  
194       TH1F *fhZPCSpectrum  = (TH1F*)l->FindObject("fhZPCSpectrum");  
195       TH1F *fhZPASpectrum  = (TH1F*)l->FindObject("fhZPASpectrum");  
196       TH1F *fhZEM1Spectrum = (TH1F*)l->FindObject("fhZEM1Spectrum"); 
197       TH1F *fhZEM2Spectrum = (TH1F*)l->FindObject("fhZEM2Spectrum"); 
198       /*TH1F *fhZNCpmc   = (TH1F*)l->FindObject("fhZNCpmc");   
199       TH1F *fhZNApmc   = (TH1F*)l->FindObject("fhZNApmc");   
200       TH1F *fhZPCpmc   = (TH1F*)l->FindObject("fhZPCpmc");   
201       TH1F *fhZPApmc   = (TH1F*)l->FindObject("fhZPApmc");*/   
202       TH2F *fhZNCCentroid  = (TH2F*)l->FindObject("fhZNCCentroid");  
203       TH2F *fhZNACentroid  = (TH2F*)l->FindObject("fhZNACentroid");   
204       TH1F *fhZNCemd   = (TH1F*)l->FindObject("fhZNCemd");   
205       TH1F *fhZNAemd   = (TH1F*)l->FindObject("fhZNAemd");   
206       TH1F *fhPMCZNCemd    = (TH1F*)l->FindObject("fhPMCZNCemd");    
207       TH1F *fhPMCZNAemd    = (TH1F*)l->FindObject("fhPMCZNAemd");    
208  
209       TH1D *hxZNC = fhZNCCentroid->ProjectionX("hxZNC");
210       TH1D *hyZNC = fhZNCCentroid->ProjectionY("hyZNC");
211       TH1D *hxZNA = fhZNACentroid->ProjectionX("hxZNA");
212       TH1D *hyZNA = fhZNACentroid->ProjectionY("hyZNA");
213       
214       Int_t index = 0;
215       xnt[index++] = (Float_t)iRun;
216       xnt[index++] = fhZNCSpectrum->GetMean();
217       xnt[index++] = fhZPCSpectrum->GetMean();
218       xnt[index++] = fhZNASpectrum->GetMean();
219       xnt[index++] = fhZPASpectrum->GetMean();
220       xnt[index++] = fhZEM1Spectrum->GetMean();
221       xnt[index++] = fhZEM2Spectrum->GetMean();
222       if(fhZNCSpectrum->GetEntries()>0)  xnt[index++] = fhZNCSpectrum->GetMean()/fhZNCSpectrum->GetEntries();
223       if(fhZPCSpectrum->GetEntries()>0)  xnt[index++] = fhZPCSpectrum->GetMean()/fhZPCSpectrum->GetEntries();
224       if(fhZNASpectrum->GetEntries()>0)  xnt[index++] = fhZNASpectrum->GetMean()/fhZNASpectrum->GetEntries();
225       if(fhZPASpectrum->GetEntries()>0)  xnt[index++] = fhZPASpectrum->GetMean()/fhZPASpectrum->GetEntries();
226       if(fhZEM1Spectrum->GetEntries()>0) xnt[index++] = fhZEM1Spectrum->GetMean()/fhZEM1Spectrum->GetEntries();
227       if(fhZEM2Spectrum->GetEntries()>0) xnt[index++] = fhZEM2Spectrum->GetMean()/fhZEM2Spectrum->GetEntries();
228       xnt[index++] = fhZNCemd->GetMean();
229       xnt[index++] = fhZNAemd->GetMean();
230       if(fhZNCemd->GetEntries()>0) xnt[index++] = fhZNCemd->GetMean()/fhZNCemd->GetEntries();
231       if(fhZNAemd->GetEntries()>0) xnt[index++] = fhZNAemd->GetMean()/fhZNAemd->GetEntries();
232       xnt[index++] = fhPMCZNCemd->GetMean();
233       xnt[index++] = fhPMCZNAemd->GetMean();
234       if(fhPMCZNCemd->GetEntries()>0) xnt[index++] = fhPMCZNCemd->GetMean()/fhPMCZNCemd->GetEntries();
235       if(fhPMCZNAemd->GetEntries()>0) xnt[index++] = fhPMCZNAemd->GetMean()/fhPMCZNAemd->GetEntries();
236       xnt[index++] = hxZNC->GetMean();
237       xnt[index++] = hyZNC->GetMean();
238       xnt[index++] = hxZNA->GetMean();
239       xnt[index++] = hyZNA->GetMean();
240       xnt[index++] = fhTDCZNSum->GetMean();
241       xnt[index++] = fhTDCZNDiff->GetMean();
242       ntzdc->Fill(xnt);
243     }
244   }
245   printf("Number of analyzed files  =  %d\n",nAnalyzedFiles);
246
247   if(nAnalyzedFiles>0){
248     TFile* outfil = new TFile(outFilNam.Data(),"recreate");
249     outfil->cd();
250     ntzdc->Write();
251     outfil->Close();
252     
253     MakePlots(outFilNam);
254   }
255 }
256
257 void MakePlots(TString ntupleFileName){
258   TFile* fil = new TFile(ntupleFileName.Data(),"read");
259   if(!fil){
260     printf("File with ntuple does not exist\n");
261     return;
262   }
263   TNtuple* ntzdc = (TNtuple*)fil->Get("ntzdc");
264
265   Float_t nrun;
266   Float_t meanZNC,meanZPC,meanZNA,meanZPA,meanZEM1,meanZEM2;
267   Float_t emeanZNC,emeanZPC,emeanZNA,emeanZPA,emeanZEM1,emeanZEM2;
268   Float_t ZNClg,ZNAlg,pmcZNClg,pmcZNAlg;
269   Float_t eZNClg,eZNAlg,epmcZNClg,epmcZNAlg;
270   Float_t xZNC,yZNC,xZNA,yZNA,tdcSum,tdcDiff;
271   
272   ntzdc->SetBranchAddress("nrun",&nrun);
273   ntzdc->SetBranchAddress("meanZNC",&meanZNC);
274   ntzdc->SetBranchAddress("meanZPC",&meanZPC);
275   ntzdc->SetBranchAddress("meanZNA",&meanZNA);
276   ntzdc->SetBranchAddress("meanZPA",&meanZPA);
277   ntzdc->SetBranchAddress("meanZEM1",&meanZEM1);
278   ntzdc->SetBranchAddress("meanZEM2",&meanZEM2);
279   ntzdc->SetBranchAddress("emeanZNC",&emeanZNC);
280   ntzdc->SetBranchAddress("emeanZPC",&emeanZPC);
281   ntzdc->SetBranchAddress("emeanZNA",&emeanZNA);
282   ntzdc->SetBranchAddress("emeanZPA",&emeanZPA);
283   ntzdc->SetBranchAddress("emeanZEM1",&emeanZEM1);
284   ntzdc->SetBranchAddress("emeanZEM2",&emeanZEM2);
285   ntzdc->SetBranchAddress("ZNClg",&ZNClg);
286   ntzdc->SetBranchAddress("ZNAlg",&ZNAlg);
287   ntzdc->SetBranchAddress("eZNClg",&eZNClg);
288   ntzdc->SetBranchAddress("eZNAlg",&eZNAlg);
289   ntzdc->SetBranchAddress("pmcZNClg",&pmcZNClg);
290   ntzdc->SetBranchAddress("pmcZNAlg",&pmcZNAlg);
291   ntzdc->SetBranchAddress("epmcZNClg",&epmcZNClg);
292   ntzdc->SetBranchAddress("epmcZNAlg",&epmcZNAlg);
293   ntzdc->SetBranchAddress("xZNC",&xZNC);
294   ntzdc->SetBranchAddress("yZNC",&yZNC);
295   ntzdc->SetBranchAddress("xZNA",&xZNA);
296   ntzdc->SetBranchAddress("yZNA",&yZNA);
297   ntzdc->SetBranchAddress("tdcSum",&tdcSum);
298   ntzdc->SetBranchAddress("tdcDiff",&tdcDiff);
299
300   TH1F *hznc      = new TH1F("hznc","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
301   TH1F *hzna      = new TH1F("hzna","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
302   TH1F *hzpc      = new TH1F("hzpc","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
303   TH1F *hzpa      = new TH1F("hzpa","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
304   TH1F *hzem1     = new TH1F("hzem1","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
305   TH1F *hzem2     = new TH1F("hzem2","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
306   TH1F *hznclg    = new TH1F("hznclg","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
307   TH1F *hznalg    = new TH1F("hznalg","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
308   TH1F *hpmcznclg = new TH1F("hpmcznclg","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
309   TH1F *hpmcznalg = new TH1F("hpmcznalg","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
310   TH1F *hxznc     = new TH1F("hxznc","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
311   TH1F *hyznc     = new TH1F("hyznc","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
312   TH1F *hxzna     = new TH1F("hxzna","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
313   TH1F *hyzna     = new TH1F("hyzna","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
314   TH1F *htdcsum   = new TH1F("htdcsum","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
315   TH1F *htdcdiff  = new TH1F("htdcdiff","",(Int_t)ntzdc->GetEntries(),0.,ntzdc->GetEntries());
316
317   for(Int_t i = 0; i<ntzdc->GetEntries();i++){
318     ntzdc->GetEvent(i);
319     //
320     hznc->SetBinContent(i+1, meanZNC);
321     hznc->SetBinError(i+1, emeanZNC);
322     hznc->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
323     hznc->GetXaxis()->SetTitle("RUN #");
324     hznc->GetYaxis()->SetTitle("ZNC mean signal");
325     hzpc->SetBinContent(i+1, meanZPC);
326     hzpc->SetBinError(i+1, emeanZPC);
327     hzpc->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
328     hzpc->GetXaxis()->SetTitle("RUN #");
329     hzpc->GetYaxis()->SetTitle("ZPC mean signal");
330     hzna->SetBinContent(i+1, meanZNA);
331     hzna->SetBinError(i+1, emeanZNA);
332     hzna->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
333     hzna->GetXaxis()->SetTitle("RUN #");
334     hzna->GetYaxis()->SetTitle("ZNA mean signal");
335     hzpa->SetBinContent(i+1, meanZPA);
336     hzpa->SetBinError(i+1, emeanZPA);
337     hzpa->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
338     hzpa->GetXaxis()->SetTitle("RUN #");
339     hzpa->GetYaxis()->SetTitle("ZPA mean signal");
340     hzem1->SetBinContent(i+1, meanZEM1);
341     hzem1->SetBinError(i+1, emeanZEM1);
342     hzem1->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
343     hzem1->GetXaxis()->SetTitle("RUN #");
344     hzem1->GetYaxis()->SetTitle("ZEM1 mean signal");
345     hzem2->SetBinContent(i+1, meanZEM2);
346     hzem2->SetBinError(i+1, emeanZEM2);
347     hzem2->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
348     hzem2->GetXaxis()->SetTitle("RUN #");
349     hzem2->GetYaxis()->SetTitle("ZEM1 mean signal");
350     hznclg->SetBinContent(i+1, ZNClg);
351     hznclg->SetBinError(i+1, eZNClg);
352     hznclg->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
353     hznclg->GetXaxis()->SetTitle("RUN #");
354     hznclg->GetYaxis()->SetTitle("ZNC LG mean signal");
355     hznalg->SetBinContent(i+1, ZNAlg);
356     hznalg->SetBinError(i+1, eZNAlg);
357     hznalg->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
358     hznalg->GetXaxis()->SetTitle("RUN #");
359     hznalg->GetYaxis()->SetTitle("ZNA LG mean signal");
360     hpmcznclg->SetBinContent(i+1, pmcZNClg);
361     hpmcznclg->SetBinError(i+1, epmcZNClg);
362     hpmcznclg->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
363     hpmcznclg->GetXaxis()->SetTitle("RUN #");
364     hpmcznclg->GetYaxis()->SetTitle("ZNC LG PMC mean signal");
365     hpmcznalg->SetBinContent(i+1, pmcZNAlg);
366     hpmcznalg->SetBinError(i+1, epmcZNAlg);
367     hpmcznalg->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
368     hpmcznalg->GetXaxis()->SetTitle("RUN #");
369     hpmcznalg->GetYaxis()->SetTitle("ZNC LG PMC mean signal");
370     hxznc->SetBinContent(i+1, xZNC);
371     hxznc->SetBinError(i+1, 0.1*xZNC);
372     hxznc->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
373     hxznc->GetXaxis()->SetTitle("RUN #");
374     hxznc->GetYaxis()->SetTitle("X_{ZN} (cm)");
375     hyznc->SetBinContent(i+1, yZNC);
376     hyznc->SetBinError(i+1, 0.1*yZNC);
377     hyznc->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
378     hyznc->GetXaxis()->SetTitle("RUN #");
379     hyznc->GetYaxis()->SetTitle("Y_{ZN} (cm)");
380     hxzna->SetBinContent(i+1, xZNA);
381     hxzna->SetBinError(i+1, 0.1*xZNA);
382     hxzna->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
383     hxzna->GetXaxis()->SetTitle("RUN #");
384     hxzna->GetYaxis()->SetTitle("X_{ZN} (cm)");
385     hyzna->SetBinContent(i+1, yZNA);
386     hyzna->SetBinError(i+1, 0.1*yZNA);
387     hyzna->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
388     hyzna->GetXaxis()->SetTitle("RUN #");
389     hyzna->GetYaxis()->SetTitle("Y_{ZN} (cm)");
390     htdcsum->SetBinContent(i+1, tdcSum);
391     htdcsum->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
392     htdcsum->GetXaxis()->SetTitle("RUN #");
393     htdcsum->GetYaxis()->SetTitle("TDC Sum (ns)");
394     htdcdiff->SetBinContent(i+1, tdcDiff);
395     htdcdiff->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
396     htdcdiff->GetXaxis()->SetTitle("RUN #");
397     htdcdiff->GetYaxis()->SetTitle("TDC Diff (ns)");
398   }
399
400   
401   TCanvas *c1 = new TCanvas("c1", "Mean value ZNs", 0, 0, 1200, 1000);
402   c1->Divide(1,2);
403   //
404   c1->cd(1); 
405   hznc->SetMarkerColor(kAzure+6); hznc->SetLineColor(kAzure+6); 
406   hznc->SetMarkerStyle(21);
407   hznc->Draw(""); 
408   hzna->SetMarkerColor(kPink-2);  hzna->SetLineColor(kPink-2);
409   hzna->SetMarkerStyle(20);
410   hzna->Draw("SAME"); 
411   //
412   TLegend *l1 = new TLegend(0.44,0.18,0.54,0.32);
413   l1->SetFillColor(kWhite);
414   l1->AddEntry(hznc," ZNC " ,"P");
415   l1->AddEntry(hzna," ZNA " ,"P");
416   l1->Draw("");
417   //
418   c1->cd(2);
419   hznalg->SetMarkerColor(kBlue+1); hznalg->SetLineColor(kBlue+1); 
420   hznalg->SetMarkerStyle(20);  hznalg->SetMinimum(0);
421   hznalg->Draw("");
422   hznclg->SetMarkerColor(kPink+5); hznclg->SetLineColor(kPink+5); 
423   hznclg->SetMarkerStyle(21);
424   hznclg->Draw("SAME"); 
425   //
426   TLegend *l2 = new TLegend(0.44,0.18,0.54,0.32);
427   l2->SetFillColor(kWhite);
428   l2->AddEntry(hznc," ZNC LG " ,"P");
429   l2->AddEntry(hzna," ZNA LG " ,"P");
430   l2->Draw("");
431   
432
433   TCanvas *c1b = new TCanvas("c1b", "Mean value ZEMs ZPs", 200, 0, 1200, 1000);
434   c1b->Divide(1,2);
435   c1b->cd(1); 
436   hzpc->SetMarkerColor(kAzure+6); hzpc->SetLineColor(kAzure+6); 
437   hzpc->SetMarkerStyle(21);
438   hzpc->Draw(""); 
439   hzpa->SetMarkerColor(kPink-2);  hzpa->SetLineColor(kPink-2);
440   hzpa->SetMarkerStyle(20);
441   hzpa->Draw("SAME"); 
442   //
443   TLegend *l3 = new TLegend(0.44,0.18,0.54,0.32);
444   l3->SetFillColor(kWhite);
445   l3->AddEntry(hzpc," ZPC " ,"P");
446   l3->AddEntry(hzpa," ZPA " ,"P");
447   l3->Draw("");
448   
449   c1b->cd(2);
450   hzem1->SetMarkerColor(kTeal-7); hzem1->SetLineColor(kTeal-7);
451   hzem1->SetMarkerStyle(29);
452   hzem2->SetMarkerColor(kTeal+5); hzem2->SetLineColor(kTeal+5);
453   hzem2->SetMarkerStyle(30);
454   hzem2->Draw(""); 
455   hzem1->Draw("SAME"); 
456   //
457   TLegend *l4 = new TLegend(0.44,0.18,0.54,0.32);
458   l4->SetFillColor(kWhite);
459   l4->AddEntry(hzem1," ZEM1 " ,"P");
460   l4->AddEntry(hzem2," ZEM2 " ,"P");
461   l4->Draw("");
462  
463   
464   /*TCanvas *c2 = new TCanvas("c2", "ZN centroids", 400, 400, 1200, 1000);
465   c2->Divide(1,2);
466   //
467   c2->cd(1);
468   hxznc->SetMarkerColor(kAzure+5); hxznc->SetLineColor(kAzure+5);
469   hxznc->SetMarkerStyle(21); 
470   hxzna->SetMarkerColor(kPink+5); hxzna->SetLineColor(kPink+5);
471   hxzna->SetMarkerStyle(20); 
472   hxznc->Draw(""); 
473   hxzna->Draw("SAME");
474   c2->cd(2);
475   hyznc->SetMarkerColor(kAzure+5); hyznc->SetLineColor(kAzure+5);
476   hyznc->SetMarkerStyle(21); 
477   hyzna->SetMarkerColor(kPink+5); hyzna->SetLineColor(kPink+5);
478   hyzna->SetMarkerStyle(20); 
479   hyznc->Draw(""); 
480   hyzna->Draw("SAME");*/
481
482   
483   TCanvas *c3 = new TCanvas("c3", "LG signals", 600, 0, 1200, 1000);
484   c3->Divide(1,2);
485   //
486   c3->cd(1);
487   hznalg->SetMarkerColor(kBlue+1); hznalg->SetLineColor(kBlue+1); 
488   hznalg->SetMarkerStyle(20);
489   hznalg->Draw(""); 
490   hznclg->SetMarkerColor(kPink+5); hznclg->SetLineColor(kPink+5); 
491   hznclg->SetMarkerStyle(21);
492   hznclg->Draw("SAME"); 
493   TLegend *l7 = new TLegend(0.44,0.18,0.54,0.32);
494   l7->SetFillColor(kWhite);
495   l7->AddEntry(hznclg," ZNC LG " ,"P");
496   l7->AddEntry(hznalg," ZNA LG " ,"P");
497   l7->Draw("");
498   //
499   c3->cd(2);
500   hpmcznalg->SetMarkerColor(kAzure+7); hpmcznalg->SetLineColor(kAzure+7); 
501   hpmcznalg->SetMarkerStyle(20); hpmcznalg->SetMinimum(0);
502   hpmcznalg->Draw(""); 
503   hpmcznclg->SetMarkerColor(kPink+6); hpmcznclg->SetLineColor(kPink+6); 
504   hpmcznclg->SetMarkerStyle(21);
505   hpmcznclg->Draw("SAME"); 
506   TLegend *l8 = new TLegend(0.44,0.15,0.58,0.32);
507   l8->SetFillColor(kWhite);
508   l8->AddEntry(hpmcznclg," ZNC LG PMC" ,"P");
509   l8->AddEntry(hpmcznalg," ZNA LG PMC" ,"P");
510   l8->Draw("");
511 }