]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/VZERO/readQAPbPb.C
Fixed end-of-line format
[u/mrichter/AliRoot.git] / PWGPP / VZERO / readQAPbPb.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2   #include <TError.h>
3   #include <TROOT.h>
4   #include <TKey.h>
5   #include <TH2.h>
6   #include <TF1.h>
7   #include <TFile.h>
8   #include <TCanvas.h>
9   #include <TPad.h>
10   #include <TStyle.h>
11   #include <TGrid.h>
12   #include <TGridResult.h>
13   #include <TEnv.h>
14   #include <TLegend.h>
15   #include <TMath.h>
16   #include <TSpectrum.h>
17
18   #include "AliCDBManager.h"
19   #include "AliCDBEntry.h"
20   #include "AliGRPObject.h"
21   #include "AliTriggerInput.h"
22   #include "AliTriggerConfiguration.h"
23 #endif
24
25 void readQAPbPb(const char * period ="LHC11h", char * files = "list.txt"){
26         gStyle->SetPalette(1);
27          TGrid::Connect("alien://");
28         // gSystem->Exec("alien_find /alice/data/2011/LHC11a/* ESDs/pass1/QA*/QAresults.root > list.txt");
29         // gSystem->Exec("sed '$d' < list.txt > tmplist.txt ; mv tmplist.txt list.txt");
30
31    // TGridResult *res = gGrid->Query(Form("/alice/data/2011/%s/*",period),"ESDs/pass1/QAresults.root");
32    // const Int_t nFiles = res->GetEntries();
33    // if (nFiles ==0) {
34    //   Error("QA","No QA files found");
35    //   delete res;
36    //   return;
37    // }
38
39         TArrayI valid(1);
40         TArrayI runs(1);
41         AliCDBManager *man = AliCDBManager::Instance();
42
43   man->SetDefaultStorage("raw://");
44
45 Int_t minRun = 1000000;
46 Int_t maxRun = 0;
47 Int_t nValid =0;
48
49 FILE * fin = fopen(files,"r");
50 Int_t runNumber, nfiles=0;
51 while(EOF!=fscanf(fin,"%d\n",&runNumber)){
52         valid.Set(nfiles+1);
53         runs.Set(nfiles+1);
54         valid.SetAt(kTRUE,nfiles);
55         runs.SetAt(runNumber,nfiles);
56
57
58  gEnv->SetValue("XNet.ConnectTimeout",10);
59  gEnv->SetValue("XNet.RequestTimeout",10);
60  gEnv->SetValue("XNet.MaxRedirectCount",2);
61  gEnv->SetValue("XNet.ReconnectTimeout",10);
62  gEnv->SetValue("XNet.FirstConnectMaxCnt",1);   
63
64         man->SetRun(runNumber);
65         AliCDBEntry *entry2=0;
66         entry2 = man->Get("GRP/GRP/Data");
67         AliGRPObject* fGRPData=0;
68         if (entry2) {
69         printf("Found an AliGRPObject in GRP/GRP/Data, reading it\n");
70         fGRPData = dynamic_cast<AliGRPObject*>(entry2->GetObject());  // new GRP entry
71         entry2->SetOwner(0);
72         }
73         TString activeDetList(AliDAQ::ListOfTriggeredDetectors(fGRPData->GetDetectorMask()));
74         TString runType(fGRPData->GetRunType());
75         TString beamType(fGRPData->GetBeamType());
76         TString machineMode(fGRPData->GetMachineMode());
77         TString lhcState(fGRPData->GetLHCState());
78         cout<<"beamType "<<beamType<<endl;
79         cout<<"machineMode "<<machineMode<<endl;
80         cout<<"lhcState "<<lhcState<<endl;
81         
82         time_t duration = fGRPData->GetTimeEnd() - fGRPData->GetTimeStart();
83
84         if(!lhcState.Contains("STABLE BEAMS")){ // Remove no BEAM runs
85                 valid.SetAt(kFALSE,nfiles);
86                 continue;
87         }
88         if(!activeDetList.Contains("VZERO")){ // Remove Runs where VZERO is not active
89                 valid.SetAt(kFALSE,nfiles);
90                 continue;
91         }
92         if(!runType.Contains("PHYSICS")){ // Remove no PHYSICS runs
93                 valid.SetAt(kFALSE,nfiles);
94                 continue;
95         }
96         if(duration<600){ // Remove Runs shorter than 10 min
97                 valid.SetAt(kFALSE,nfiles);
98                 continue;
99         }
100         nfiles++;
101         nValid++;
102         
103         if(runNumber>maxRun) maxRun = runNumber;
104         if(runNumber<minRun) minRun = runNumber;
105 }
106
107 TH1F * hTimeA = new TH1F("hTimeA","BB Leading time;;Time (ns)",nValid,-0.5,nValid-0.5);
108 TH1F * hTimeC = new TH1F("hTimeC","BB Leading time;;Time (ns)",nValid,-0.5,nValid-0.5);
109 TH1F * hBB_BG = new TH1F("hBB_BG","Trigger ratio",nValid,-0.5,nValid-0.5);
110 TH1F * hBB_EE = new TH1F("hBB_EE","Trigger ratio",nValid,-0.5,nValid-0.5);
111 TH1F * hAdcA = new TH1F("hAdcA","Average Charge",nValid,-0.5,nValid-0.5);
112 TH1F * hAdcC = new TH1F("hAdcC","Average Charge",nValid,-0.5,nValid-0.5);
113 TH1F * hMultA = new TH1F("hMultA","Average number of Fired cell",nValid,-0.5,nValid-0.5);
114 TH1F * hMultC = new TH1F("hMultC","Average number of Fired cell",nValid,-0.5,nValid-0.5);
115 TH1F * hTriggerEff_CVLN = new TH1F("hTriggerEff_CVLN","CVLN / CVBN",nValid,-0.5,nValid-0.5);
116 TH1F * hTriggerEff_CVHN = new TH1F("hTriggerEff_CVHN","CVHN / CVBN",nValid,-0.5,nValid-0.5);
117 TH1F * hTriggerEff_CVHN2 = new TH1F("hTriggerEff_CVHN2","CVHN / CVLN",nValid,-0.5,nValid-0.5);
118 TH1F * hPMTEdges[64];
119 for(int i = 0; i < 64; ++i){
120         hPMTEdges[i] = new  TH1F(Form("hPMTEdges%d",i),Form("Multiplicity edge Cell %d",i),nValid,-0.5,nValid-0.5);
121 }
122
123 int nEntries=0;
124 TString trigMB, trigCVLN, trigCVHN;
125 for(int ifile = nfiles-1; ifile >= 0; ifile--){
126
127         if(!valid.At(ifile)) continue;
128         runNumber = runs.At(ifile);
129
130         if(runNumber>=169683) {
131                 trigMB   = "CPBI2_B1";
132                 trigCVLN = "CSEMI_R1";
133                 trigCVHN = "CCENT_R2";
134         } else if(runNumber>168171){
135                 trigMB   = "CPBI2_B1";
136                 trigCVLN = "CVLN_R1";
137                 trigCVHN = "CVHN_R2";
138         } else if(runNumber>167693){
139                 trigMB   = "CPBI2_B1";
140                 trigCVLN = "CVLN_B2";
141                 trigCVHN = "CVHN_R2";
142         } else if(runNumber>166532)     {
143                 trigMB   = "CPBI1";
144                 trigCVLN = "CVLN";
145                 trigCVHN = "CVHN";
146         } else {
147                 trigMB   = "CPBI1";
148                 trigCVLN = "CVLN";
149                 trigCVHN = "CVHN";
150         }
151                 
152         TGridResult *res;
153         if(strcmp(period,"LHC11h")==0){
154                 res = gGrid->Query(Form("/alice/data/2011/%s/000%d/*",period,runNumber),"ESDs/pass1_HLT/QAresults.root");
155         }
156
157
158 if (res->GetEntries() ==0) {
159      Error("QA",Form("No QA files found for run %d\n",runNumber));
160      delete res;
161          nEntries++;
162      continue;
163    }
164
165            man->SetRun(runNumber);
166         
167         AliCDBEntry *entryCTP = man->Get("GRP/CTP/Config");
168    AliTriggerConfiguration *configCTP = (AliTriggerConfiguration*)entryCTP->GetObject();
169         TObjArray  inputsArray = configCTP->GetInputs();
170         
171         Double_t rnd1=1., rnd2=1., bc1=1., bc2=1.;
172         AliTriggerInput * input;
173         input = (AliTriggerInput*)(inputsArray.FindObject("RND1"));
174         if(input) rnd1 =  (input->GetSignature())/(double)(0x7fffffff );
175         if(input)cout<<Form("RND1 = %d",input->GetSignature())<<endl;
176
177         input = (AliTriggerInput*)(inputsArray.FindObject("RND2"));
178         if(input) rnd2 =  (input->GetSignature())/(double)(0x7fffffff );
179         if(input)cout<<Form("RND2 = %d",input->GetSignature())<<endl;
180         
181         input = (AliTriggerInput*)(inputsArray.FindObject("BC1"));
182         if(input) bc1 =  1./(input->GetSignature()+1.);
183         if(input)cout<<Form("BC1 = %d",input->GetSignature())<<endl;
184         
185         input = (AliTriggerInput*)(inputsArray.FindObject("BC2"));
186         if(input) bc2 =  1./(input->GetSignature()+1.);
187         if(input)cout<<Form("BC2 = %d",input->GetSignature())<<endl;
188
189
190         
191         TString filename = res->GetKey(0, "turl");
192         //      TObjArray* tmp = filename.Tokenize("/");
193         if(filename == "") continue;
194         TFile *fQA = TFile::Open(filename.Data());
195         if (!fQA) {
196         Error("QA",Form("Can not open QA file found for run %d\n",runNumber));
197                 nEntries++;
198         continue;
199         }
200         TList *list = (TList*)fQA->Get("VZERO_PbPb_Performance/PbPbVZEROHists");
201         TH2F *hTriggerDecision = (TH2F*)list->FindObject("hTriggerDecision");
202         TH1F *hAdcNoTimeA = (TH1F*)list->FindObject("hAdcNoTimeV0A");
203         TH1F *hAdcWithTimeA = (TH1F*)list->FindObject("hAdcWithTimeV0A");
204         TH1F *hAdcNoTimeC = (TH1F*)list->FindObject("hAdcNoTimeV0C");
205         TH1F *hAdcWithTimeC = (TH1F*)list->FindObject("hAdcWithTimeV0C");
206         TH2F *hadcpmtwithtime = (TH2F*)list->FindObject("hadcpmtwithtime");     
207         TH1F *htimepmtA = (TH1F*)list->FindObject("htimepmtV0A");
208         TH1F *htimepmtC = (TH1F*)list->FindObject("htimepmtV0C");
209         TH1F *hwidthA = (TH1F*)list->FindObject("hwidthV0A");
210         TH1F *hwidthC = (TH1F*)list->FindObject("hwidthV0C");
211         //      TH1F *hV0ampl = (TH1F*)list->FindObject("hV0ampl");
212         TH2F *htimepmt = (TH2F*)list->FindObject("htimepmt");   
213         TH2F *hwidthpmt = (TH2F*)list->FindObject("hwidthpmt"); 
214         TH2F *hadcwidthA = (TH2F*)list->FindObject("hadcwidthV0A");     
215         TH2F *hadcwidthC = (TH2F*)list->FindObject("hadcwidthV0C");     
216         TH2F *hAdcTimeA = (TH2F*)list->FindObject("hAdcTimeV0A");       
217         TH2F *hAdcTimeC = (TH2F*)list->FindObject("hAdcTimeV0C");       
218         TH2F *htimecorr = (TH2F*)list->FindObject("htimecorr"); 
219         TH2F *hNFlags   = (TH2F*)list->FindObject("hNFlags");
220         TH1F *hV0A = (TH1F*)hNFlags->ProjectionX("hV0A",1,hNFlags->GetNbinsY());
221         TH1F *hV0C = (TH1F*)hNFlags->ProjectionY("hV0C",1,hNFlags->GetNbinsX());
222         TH2F* hVtxXYBB  =(TH2F*) list->FindObject("fhVtxXYBB");
223         TH1F* hVtxZBB   =(TH1F*) list->FindObject("fhVtxZBB");
224         TH2F* hVtxXYBGA =(TH2F*) list->FindObject("fhVtxXYBGA");
225         TH1F* hVtxZBGA  =(TH1F*) list->FindObject("fhVtxZBGA");
226         TH2F* hVtxXYBGC =(TH2F*) list->FindObject("fhVtxXYBGC");
227         TH1F* hVtxZBGC  =(TH1F*) list->FindObject("fhVtxZBGC");
228         
229         TH2F* hRecoMult = (TH2F*) list->FindObject(Form("hRecoMult_%s-",trigMB.Data()));
230         TH2F* hRecoMultPMT = (TH2F*) list->FindObject(Form("hRecoMultPMT_%s-",trigMB.Data()));
231         TH1F* hTotRecoMult = (TH1F*) list->FindObject(Form("hTotRecoMult_%s-",trigMB.Data()));
232         TH1F* hTotRecoMult_CVLN = (TH1F*) list->FindObject(Form("hTotRecoMult_%s-",trigCVLN.Data()));
233         TH1F* hTotRecoMult_CVHN = (TH1F*) list->FindObject(Form("hTotRecoMult_%s-",trigCVHN.Data()));
234         TH2F* hEqualizedMult = (TH2F*) list->FindObject(Form("hEqualizedMult_%s-",trigMB.Data()));
235         
236         Double_t BB  = hTriggerDecision->GetBinContent(2,2);
237         Double_t EE  = hTriggerDecision->GetBinContent(1,1);
238         Double_t BGA = hTriggerDecision->GetBinContent(3,2);
239         Double_t BGC = hTriggerDecision->GetBinContent(2,3);
240                 
241         hTimeA->GetXaxis()->SetBinLabel(nEntries+1,Form("%d",runNumber));
242         hTimeC->GetXaxis()->SetBinLabel(nEntries+1,Form("%d",runNumber));
243         hBB_BG->GetXaxis()->SetBinLabel(nEntries+1,Form("%d",runNumber));
244         hBB_EE->GetXaxis()->SetBinLabel(nEntries+1,Form("%d",runNumber));
245         hAdcA->GetXaxis()->SetBinLabel(nEntries+1,Form("%d",runNumber));
246         hAdcC->GetXaxis()->SetBinLabel(nEntries+1,Form("%d",runNumber));
247         hMultA->GetXaxis()->SetBinLabel(nEntries+1,Form("%d",runNumber));
248         hMultC->GetXaxis()->SetBinLabel(nEntries+1,Form("%d",runNumber));
249         hTriggerEff_CVLN->GetXaxis()->SetBinLabel(nEntries+1,Form("%d",runNumber));
250         hTriggerEff_CVHN->GetXaxis()->SetBinLabel(nEntries+1,Form("%d",runNumber));
251         hTriggerEff_CVHN2->GetXaxis()->SetBinLabel(nEntries+1,Form("%d",runNumber));
252     for(int i = 0; i < 64; ++i) {
253         hPMTEdges[i]->GetXaxis()->SetBinLabel(nEntries+1,Form("%d",runNumber));
254     }
255         
256         if(hAdcWithTimeA->GetEntries()==0) {
257                 delete fQA;     
258                 nEntries++;             
259                 continue;
260         }
261
262         Double_t cVLN = hTotRecoMult_CVLN->GetEntries();
263         Double_t cVHN = hTotRecoMult_CVHN->GetEntries();
264         Double_t cVBN = hTotRecoMult->GetEntries();
265         
266         // Correct for downscaling factor
267         Double_t scaleVBN = 1., scaleVLN = 1., scaleVHN = 1.;
268         if(trigMB.Contains("B1")) {
269                 if(bc1>0.) scaleVBN=1./bc1;
270         } else  if(trigMB.Contains("B2")) {
271                 if(bc2>0.) scaleVBN=1./bc2;
272         } else if(trigMB.Contains("R1")) {
273                 if(rnd1>0.) scaleVBN=1./rnd1;
274         } else if(trigMB.Contains("R2")) {
275                 if(rnd2>0.) scaleVBN=1./rnd2;
276         }
277         if(trigCVLN.Contains("B1")) {
278                 if(bc1>0.) scaleVLN=1./bc1;
279         } else  if(trigCVLN.Contains("B2")) {
280                 if(bc2>0.) scaleVLN=1./bc2;
281         } else if(trigCVLN.Contains("R1")) {
282                 if(rnd1>0.) scaleVLN=1./rnd1;
283         } else if(trigCVLN.Contains("R2")) {
284                 if(rnd2>0.) scaleVLN=1./rnd2;
285         }
286         
287         if(trigCVHN.Contains("B1")) {
288                 if(bc1>0.) scaleVHN=1./bc1;
289         } else  if(trigCVHN.Contains("B2")) {
290                 if(bc2>0.) scaleVHN=1./bc2;
291         } else if(trigCVHN.Contains("R1")) {
292                 if(rnd1>0.) scaleVHN=1./rnd1;
293         } else if(trigCVHN.Contains("R2")) {
294                 if(rnd2>0.) scaleVHN=1./rnd2;
295         }
296
297         cout<<Form("CVBN = %lf \nCVLN = %lf \nCVHN = %lf\n",cVBN,cVLN,cVHN);
298         hTotRecoMult->Scale(scaleVBN);
299         hTotRecoMult_CVLN->Scale(scaleVLN);
300         hTotRecoMult_CVHN->Scale(scaleVHN);
301
302         cVBN *= scaleVBN/100.;
303         cVLN *= scaleVLN;
304         cVHN *= scaleVHN;
305         
306         cout<<Form("CVBN = %lf \nCVLN = %lf \nCVHN = %lf\n",cVBN,cVLN,cVHN);
307
308         if(cVBN >0.){
309                 hTriggerEff_CVLN->SetBinContent(nEntries+1,cVBN/cVLN);
310                 if(cVLN>0.) hTriggerEff_CVLN->SetBinError(nEntries+1,cVBN/cVLN*(TMath::Sqrt(1./cVLN+1./cVBN)));
311                 hTriggerEff_CVHN->SetBinContent(nEntries+1,cVBN/cVHN);
312                 if(cVHN>0.) hTriggerEff_CVHN->SetBinError(nEntries+1,cVBN/cVHN*(TMath::Sqrt(1./(cVHN/scaleVHN)+1./(cVBN/scaleVBN*100.))));
313                 if(cVLN>0.) hTriggerEff_CVHN2->SetBinContent(nEntries+1,cVLN/cVHN);
314                 if(cVHN>0. && cVLN>0.) hTriggerEff_CVHN2->SetBinError(nEntries+1,cVLN/cVHN*(TMath::Sqrt(1./(cVHN/scaleVHN)+1./(cVLN/scaleVLN))));
315         }
316
317     Double_t beta1[64], beta2[64];
318     Double_t q = 1. - 1.e-4;
319     Double_t q2 = 1. - 2.e-4;
320     for(int i = 0; i < 64; ++i) {
321         ((TH1D*)hRecoMultPMT->ProjectionY(Form("hRecoMultPMT%d",i),i+1,i+1))->GetQuantiles(1,&beta1[i],&q);
322         ((TH1D*)hRecoMultPMT->ProjectionY(Form("hRecoMultPMT%d",i),i+1,i+1))->GetQuantiles(1,&beta2[i],&q2);
323         hPMTEdges[i]->SetBinContent(nEntries+1,(beta1[i]+beta2[i])/2.);
324                 hPMTEdges[i]->SetBinError(nEntries+1,(beta1[i] - beta2[i]));
325     }
326
327     Double_t betaSide1[2];
328     Double_t betaSide2[2];
329     for(int i = 0; i < 2; ++i) {
330                 if(i==0){
331                         ((TH1D*)hRecoMult->ProjectionY(Form("hRecoMult1%d",i)))->GetQuantiles(1,&betaSide1[i],&q);
332                         ((TH1D*)hRecoMult->ProjectionY(Form("hRecoMult2%d",i)))->GetQuantiles(1,&betaSide2[i],&q2);
333             }else{
334                         ((TH1D*)hRecoMult->ProjectionX(Form("hRecoMult1%d",i)))->GetQuantiles(1,&betaSide1[i],&q);
335                         ((TH1D*)hRecoMult->ProjectionX(Form("hRecoMult2%d",i)))->GetQuantiles(1,&betaSide2[i],&q2);
336                 }
337     }
338         hAdcA->SetBinContent(nEntries+1,(betaSide1[1] + betaSide2[1])/2.);
339         hAdcC->SetBinContent(nEntries+1,(betaSide1[0] + betaSide2[0])/2.);
340         hAdcA->SetBinError(nEntries+1,(betaSide1[1] - betaSide2[1]));
341         hAdcC->SetBinError(nEntries+1,(betaSide1[0] - betaSide2[0]));
342         
343         
344
345         
346                 TSpectrum s;
347                 int nPeaksFound;
348                 float * peaks;
349                 float max;
350                 float shiftA = 8.;
351
352                 nPeaksFound = s.Search(htimepmtA); peaks = s.GetPositionX(); max = -25.;
353                 for(int i=0;i<nPeaksFound;i++) if(peaks[i]>max) max = peaks[i]; 
354                 htimepmtA->Fit("gaus","","",max-1.,max+1.);
355                 hTimeA->Fill(nEntries,htimepmtA->GetFunction("gaus")->GetParameter(1)-shiftA);
356
357                 nPeaksFound = s.Search(htimepmtC); peaks = s.GetPositionX(); max = -25.;
358                 for(int i=0;i<nPeaksFound;i++) if(peaks[i]>max) max = peaks[i]; 
359                 htimepmtC->Fit("gaus","","",max-1.,max+1.);
360                 hTimeC->Fill(nEntries,htimepmtC->GetFunction("gaus")->GetParameter(1));
361
362                 if(BB) {
363                         hBB_BG->SetBinContent(nEntries,(BGA+BGC)/BB);
364                         hBB_EE->SetBinContent(nEntries,EE/BB);
365                 } else {
366                         hBB_BG->SetBinContent(nEntries,0);
367                         hBB_EE->SetBinContent(nEntries,0);
368                 }
369
370                 hMultA->SetBinContent(nEntries,hV0A->GetMean());
371                 hMultC->SetBinContent(nEntries,hV0C->GetMean());
372
373                 //              hAdcA->SetBinContent(nEntries,hAdcWithTimeA->GetMean());
374                 //              hAdcC->SetBinContent(nEntries,hAdcWithTimeC->GetMean());
375
376 //-------------
377         TCanvas * cOut = new TCanvas("cOut",Form("Run %d",runNumber));
378         cOut->Divide(2,2);
379         cOut->cd(1); cOut->GetPad(1)->SetLogy();
380         hAdcNoTimeA->Draw("l");
381         hAdcWithTimeA->Draw("same"); hAdcWithTimeA->SetLineColor(2);
382         
383         cOut->cd(2); cOut->GetPad(2)->SetLogy();
384         hAdcNoTimeC->Draw("l");
385         hAdcWithTimeC->Draw("same"); hAdcWithTimeC->SetLineColor(2);
386         
387         cOut->cd(3); cOut->GetPad(3)->SetLogz();
388         hadcpmtwithtime->Draw("colz");
389
390
391         cOut->cd(4); cOut->GetPad(4)->SetLogy(0);cOut->GetPad(4)->SetLogz(1);
392         hEqualizedMult->Draw("colz");  
393         
394         cOut->Print(Form("QA_Run_%d.pdf(",runNumber));
395 //-------------
396
397         cOut->Clear();
398         cOut->Divide(2,2);
399         cOut->cd(1); cOut->GetPad(1)->SetLogy();
400         htimepmtA->GetXaxis()->SetRangeUser(-25.,25.); htimepmtA->Draw();
401
402         cOut->cd(2); cOut->GetPad(2)->SetLogy();
403         htimepmtC->GetXaxis()->SetRangeUser(-25.,25.); htimepmtC->Draw();
404         
405         cOut->cd(3); cOut->GetPad(3)->SetLogy();cOut->GetPad(3)->SetLogz(0);
406         //hwidthA->GetXaxis()->SetRangeUser(0.,50.); 
407         hwidthA->Draw();
408         
409         cOut->cd(4); cOut->GetPad(4)->SetLogy();cOut->GetPad(4)->SetLogz(0);
410         //hwidthC->GetXaxis()->SetRangeUser(0.,50.); 
411         hwidthC->Draw();
412
413         cOut->Print(Form("QA_Run_%d.pdf",runNumber));
414 //-------------
415
416         cOut->Clear();
417         cOut->Divide(2,2);
418         cOut->cd(1); cOut->GetPad(1)->SetLogy(0);cOut->GetPad(1)->SetLogz();
419         htimepmt->Draw("colz");
420
421         cOut->cd(2); cOut->GetPad(2)->SetLogy(0);cOut->GetPad(2)->SetLogz();
422         //hwidthpmt->GetYaxis()->SetRangeUser(0.,50.); 
423         hwidthpmt->Draw("colz");
424         
425         cOut->cd(3); cOut->GetPad(3)->SetLogy(0);cOut->GetPad(3)->SetLogz();
426         //hadcwidthA->GetYaxis()->SetRangeUser(0.,50.); 
427         hadcwidthA->Draw("colz");
428         
429         cOut->cd(4); cOut->GetPad(4)->SetLogy(0);cOut->GetPad(4)->SetLogz();
430         //hadcwidthC->GetYaxis()->SetRangeUser(0.,50.); 
431         hadcwidthC->Draw("colz");
432
433         cOut->Print(Form("QA_Run_%d.pdf",runNumber));
434 //-------------
435
436         cOut->Clear();
437         cOut->Divide(2,2);
438         cOut->cd(1); cOut->GetPad(1)->SetLogy(0);cOut->GetPad(1)->SetLogz();
439         hAdcTimeA->Draw("colz");
440
441         cOut->cd(2); cOut->GetPad(2)->SetLogy(0);cOut->GetPad(2)->SetLogz();
442         hAdcTimeC->Draw("colz");
443         
444         cOut->cd(3); cOut->GetPad(3)->SetLogz();
445         hTriggerDecision->Draw("text");
446         
447         cOut->cd(4); cOut->GetPad(4)->SetLogz(); cOut->GetPad(4)->SetLogz(0);
448         htimecorr->Draw("colz");
449
450         cOut->Print(Form("QA_Run_%d.pdf",runNumber));
451 //-------------
452
453         cOut->Clear();
454         cOut->Divide(2,2);
455         cOut->cd(1);  cOut->GetPad(1)->SetLogy(1);cOut->GetPad(1)->SetLogz(0);
456         hV0A->GetXaxis()->SetRangeUser(0.,33.);hV0A->Draw();
457
458         cOut->cd(2); cOut->GetPad(2)->SetLogy(1);cOut->GetPad(2)->SetLogz(0);
459         hV0C->GetXaxis()->SetRangeUser(0.,33.);hV0C->Draw();
460
461         cOut->cd(3); cOut->GetPad(3)->SetLogy(0);cOut->GetPad(3)->SetLogz(0); cOut->GetPad(3)->SetGridy(1);
462         hTotRecoMult->Sumw2();
463         hTotRecoMult_CVLN->Sumw2();
464         hTotRecoMult_CVHN->Sumw2();
465         hTotRecoMult_CVLN->Divide(hTotRecoMult);
466         hTotRecoMult_CVHN->Divide(hTotRecoMult);
467         hTotRecoMult_CVLN->Draw("e"); hTotRecoMult_CVLN->SetLineColor(4);hTotRecoMult_CVLN->SetMarkerStyle(0);hTotRecoMult_CVLN->SetMaximum(1.2);
468         hTotRecoMult_CVLN->SetTitle("Multiplicity efficiency CVLN (blue) and CHVN (red)");
469         hTotRecoMult_CVHN->Draw("esame"); hTotRecoMult_CVHN->SetLineColor(2);hTotRecoMult_CVHN->SetMarkerStyle(0);
470         
471         cOut->cd(4); cOut->GetPad(4)->SetLogy(0);cOut->GetPad(4)->SetLogz(1); cOut->GetPad(4)->SetGridx(1); cOut->GetPad(4)->SetGridy(1);
472         hRecoMult->Draw("colz");
473         
474         
475         cOut->Print(Form("QA_Run_%d.pdf",runNumber));
476 //-------------
477
478         cOut->Clear();
479         cOut->Divide(2,3);
480
481         cOut->cd(1);  cOut->GetPad(1)->SetLogy(0);cOut->GetPad(1)->SetLogz(1);
482         hVtxXYBB->Draw("colz");
483
484         cOut->cd(2); cOut->GetPad(2)->SetLogy(1);cOut->GetPad(2)->SetLogz(0);
485         hVtxZBB->Draw();
486         
487         cOut->cd(3); cOut->GetPad(3)->SetLogy(0); cOut->GetPad(3)->SetLogz(1);
488         hVtxXYBGA->Draw("colz");
489         
490         cOut->cd(4); cOut->GetPad(4)->SetLogy(); cOut->GetPad(3)->SetLogz(0);
491         hVtxZBGA->Draw();
492         
493         cOut->cd(5); cOut->GetPad(5)->SetLogy(0); cOut->GetPad(5)->SetLogz(1);
494         hVtxXYBGC->Draw("colz");
495         
496         cOut->cd(6); cOut->GetPad(6)->SetLogy(); cOut->GetPad(6)->SetLogz(0);
497         hVtxZBGC->Draw();
498
499         cOut->Print(Form("QA_Run_%d.pdf)",runNumber));
500 //-------------
501
502                         delete cOut;
503                         delete fQA;                     
504
505                         nEntries++;
506
507                 }
508 //              delete res;
509
510 gStyle->SetOptStat(0);
511 hTimeA->SetMarkerStyle(20);
512 hTimeA->SetMarkerColor(2);
513
514 hTimeC->SetMarkerStyle(20);
515 hTimeC->SetMarkerColor(4);
516
517 TCanvas * c = new TCanvas("c","Leading time versus run number");
518 c->SetGridy();
519 hTimeA->Draw("P");
520 hTimeA->SetMinimum(TMath::Min(hTimeA->GetMinimum(),hTimeC->GetMinimum())-1.);
521 hTimeA->SetMaximum(TMath::Max(hTimeA->GetMaximum(),hTimeC->GetMaximum())+1.);
522 // hTimeA->GetXaxis()->SetLabelOptions("v");
523 // hTimeC->GetXaxis()->SetLabelOptions("v");
524
525 hTimeC->Draw("Psame");
526 TLegend * lg = new TLegend(0.8,0.9,1,1);
527 lg->AddEntry(hTimeA,"V0A - 8 ns","p");
528 lg->AddEntry(hTimeC,"V0C","p");
529 lg->Draw("same");
530 // TPave * pavA = new TPave(-0.5,TMath::Max(hTimeA->GetMinimum(),1.5-shiftA),nValid-0.5,TMath::Min(hTimeA->GetMaximum(),33.5-shiftA),0);
531 // pavA->SetFillStyle(3004);
532 // pavA->SetFillColor(2);
533 // TPave * pavC = new TPave(-0.5,TMath::Max(hTimeA->GetMinimum(),0.5),nValid-0.5,TMath::Min(hTimeA->GetMaximum(),25.5),0);
534 // pavC->SetFillStyle(3005);
535 // pavC->SetFillColor(4);
536 // 
537 // pavA->Draw("same");
538 // pavC->Draw("same");
539 TFile * fout = TFile::Open(Form("QA_Resume_%d_%d.root",minRun,maxRun),"RECREATE");
540
541 c->Print(Form("QA_Resume_%d_%d.pdf(",minRun,maxRun));
542 c->Write();
543
544 TCanvas * c2 = new TCanvas("c2","Trigger ratios");
545 c2->SetGridy();
546
547 hBB_BG->SetMarkerStyle(20);
548 hBB_BG->SetMarkerColor(2);
549 hBB_EE->SetMarkerStyle(20);
550 hBB_EE->SetMarkerColor(4);
551
552 hBB_BG->Draw("P");hBB_BG->SetTitle("Beam-Gas / Beam-Beam");
553 // hBB_EE->Draw("Psame");
554 // TLegend * lg2 = new TLegend(0.8,0.9,1,1);
555 // lg2->AddEntry(hBB_BG,"BG / BB","p");
556 // lg2->AddEntry(hBB_EE,"EE / BB","p");
557 // lg2->Draw("same");
558
559 c2->Print(Form("QA_Resume_%d_%d.pdf",minRun,maxRun));
560 c2->Write();
561
562
563 TCanvas * c3 = new TCanvas("c3","Multiplicity Eddes V0A/V0C");
564 c3->SetGridy();
565
566 hAdcA->SetMarkerStyle(20);
567 hAdcA->SetMarkerColor(2);
568 hAdcC->SetMarkerStyle(20);
569 hAdcC->SetMarkerColor(4);
570 hAdcA->SetMinimum(0);
571 hAdcA->SetMaximum(100);
572
573 hAdcA->Draw("P");
574 hAdcC->Draw("Psame");
575 TLegend * lg3 = new TLegend(0.8,0.9,1,1);
576 lg3->AddEntry(hAdcA,"V0A","p");
577 lg3->AddEntry(hAdcC,"V0C","p");
578 lg3->Draw("same");
579
580 c3->Print(Form("QA_Resume_%d_%d.pdf",minRun,maxRun));
581 c3->Write();
582
583
584 TCanvas * c4 = new TCanvas("c4","Average number of cell");
585 c4->SetGridy();
586
587 hMultA->SetMarkerStyle(20);
588 hMultA->SetMarkerColor(2);
589 hMultC->SetMarkerStyle(20);
590 hMultC->SetMarkerColor(4);
591 hMultA->SetMinimum(0);
592 hMultA->SetMaximum(32);
593
594 hMultA->Draw("P");
595 hMultC->Draw("Psame");
596 TLegend * lg4 = new TLegend(0.8,0.9,1,1);
597 lg4->AddEntry(hMultA,"V0A","p");
598 lg4->AddEntry(hMultC,"V0C","p");
599 lg4->Draw("same");
600
601 c4->Print(Form("QA_Resume_%d_%d.pdf",minRun,maxRun));
602 c4->Write();
603
604 TCanvas * c5 = new TCanvas("c5","Trigger Efficiency");
605 c5->SetGridy();
606
607 hTriggerEff_CVLN->SetMarkerStyle(20);
608 hTriggerEff_CVLN->SetMarkerColor(2);
609 hTriggerEff_CVHN->SetMarkerStyle(20);
610 hTriggerEff_CVHN->SetMarkerColor(4);
611 hTriggerEff_CVHN2->SetMarkerStyle(4);
612 hTriggerEff_CVHN2->SetMarkerColor(1);
613 hTriggerEff_CVLN->SetMinimum(0);
614 hTriggerEff_CVLN->SetMaximum(15.);
615 hTriggerEff_CVLN->SetTitle("Centrality triggers fraction");
616
617 hTriggerEff_CVLN->Draw("P");
618 hTriggerEff_CVHN->Draw("Psame");
619 hTriggerEff_CVHN2->Draw("Psame");
620 TLegend * lg5 = new TLegend(0.7,0.8,1,1);
621 lg5->AddEntry(hTriggerEff_CVLN,"(CPBI2/100) / CVLN","p");
622 lg5->AddEntry(hTriggerEff_CVHN,"(CPBI2/100) / CVHN","p");
623 lg5->AddEntry(hTriggerEff_CVHN2,"CVLN / CVHN","p");
624 lg5->Draw("same");
625
626 c5->Print(Form("QA_Resume_%d_%d.pdf",minRun,maxRun));
627 c5->Write();
628
629
630 TCanvas * cedge[8];
631 for(int i = 0; i < 8; ++i){
632         cedge[i] = new TCanvas(Form("cedge%d",i),Form("Edge Ring %d",i));
633         cedge[i]->SetGridy();
634         cedge[i]->Divide(3,3);
635         for(int iCh = 0; iCh < 8; ++iCh){
636                 cedge[i]->cd(iCh+1);
637                 hPMTEdges[iCh+i*8]->SetMarkerStyle(20);
638                 hPMTEdges[iCh+i*8]->Draw("P");
639         }
640         if(i==7) cedge[i]->Print(Form("QA_Resume_%d_%d.pdf)",minRun,maxRun));
641         else cedge[i]->Print(Form("QA_Resume_%d_%d.pdf",minRun,maxRun));
642         cedge[i]->Write();
643 }
644
645
646
647
648
649 fout->Close();
650
651 gSystem->Exec(Form("tar cvf QA_Runs_%d_%d.tar QA_Run*.pdf QA_Resume_%d_%d.pdf QA_Resume_%d_%d.root",minRun,maxRun,minRun,maxRun,minRun,maxRun));
652 gSystem->Exec(Form("rm -f QA_Run*.pdf QA_Resume_%d_%d.pdf  QA_Resume_%d_%d.root",minRun,maxRun,minRun,maxRun));
653
654
655 }
656