]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/macrosSDD/TrendQAtrainSDD.C
New feature to plot raw data of single module
[u/mrichter/AliRoot.git] / ITS / macrosSDD / TrendQAtrainSDD.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TH1F.h>
3 #include <TH2F.h>
4 #include <TF1.h>
5 #include <TPad.h>
6 #include <TGraphErrors.h>
7 #include <TROOT.h>
8 #include <TFile.h>
9 #include <TTree.h>
10 #include <TGrid.h>
11 #include <TGridResult.h>
12 #include <TMath.h>
13 #include <TSystem.h>
14 #include <TNtuple.h>
15 #include <TCanvas.h>
16 #include <TStyle.h>
17 #include <TLatex.h>
18 #include <TLegend.h>
19 #include <TLegendEntry.h>
20 #include "AliITSgeomTGeo.h"
21 #endif
22
23 Double_t LangausFun(Double_t *x, Double_t *par);
24 void MakePlots(TString ntupleFileName);
25
26
27 void TrendQAtrainSDD(TString period,
28                      TString recoPass,
29                      TString qaTrain="QA",
30                      Int_t firstRun=0,
31                      Int_t lastRun=999999999,
32                      TString fileName="QAresults.root"){
33
34   gStyle->SetOptStat(0);
35   
36
37   TGrid::Connect("alien:");
38   Int_t year=0;
39   if(period.Contains("LHC09")) year=2009;
40   else if(period.Contains("LHC10")) year=2010;
41   else if(period.Contains("LHC11")) year=2011;
42
43   TString outFilNam=Form("TrendingSDD_%s_%s_%s.root",period.Data(),recoPass.Data(),qaTrain.Data());
44
45
46   const Int_t nVariables=31;
47   TNtuple* ntsdd=new TNtuple("ntsdd","SDD trending","nrun:fracTrackWithClu1:errfracTrackWithClu1:fracTrackWithClu2:errfracTrackWithClu2:fracTrackWithClu3:errfracTrackWithClu3:fracTrackWithClu4:errfracTrackWithClu4:fracTrackWithClu5:errfracTrackWithClu5:fracTrackWithClu6:errfracTrackWithClu6:meanTrPts3:errmeanTrPts3:meanTrPts4:errmeanTrPts4:minDrTime:errminDrTime:meanDrTime:errmeanDrTime:fracExtra:errfracExtra:meandEdxTB0:errmeandEdxTB0:meandEdxTB5:errmeandEdxTB5:nMod95:nMod80:nMod60:nModEmpty");
48   Float_t xnt[nVariables];
49
50   TBits* readRun=new TBits(999999);
51   readRun->ResetAllBits();
52   if(!gSystem->Exec(Form("ls -l %s > /dev/null 2>&1",outFilNam.Data()))){
53     TFile* oldfil=new TFile(outFilNam.Data());
54     TNtuple* ntmp=(TNtuple*)oldfil->Get("ntsdd");
55     Bool_t isOK=kFALSE;
56     if(ntmp){
57       if(ntmp->GetNvar()==ntsdd->GetNvar()){
58         isOK=kTRUE;
59         TObjArray* arr1=(TObjArray*)ntsdd->GetListOfBranches();
60         TObjArray* arr2=(TObjArray*)ntmp->GetListOfBranches();
61         for(Int_t iV=0; iV<ntmp->GetNvar(); iV++){
62           TString vnam1=arr1->At(iV)->GetName();
63           TString vnam2=arr2->At(iV)->GetName();
64           if(vnam1!=vnam2) isOK=kFALSE;
65           ntmp->SetBranchAddress(vnam2.Data(),&xnt[iV]);
66         }
67         if(isOK){
68           for(Int_t nE=0; nE<ntmp->GetEntries(); nE++){
69             ntmp->GetEvent(nE);
70             Int_t theRun=(Int_t)(xnt[0]+0.0001);
71             readRun->SetBitNumber(theRun);
72             ntsdd->Fill(xnt);
73           }
74         }
75       }
76     }
77     if(!isOK){
78       printf("Ntuple in local file not OK -> will be recreated\n");
79     }
80     oldfil->Close();
81     delete oldfil;
82   }
83
84   if(!gGrid||!gGrid->IsConnected()) {
85     printf("gGrid not found! exit macro\n");
86     return;
87   }
88
89   TString  path=Form("/alice/data/%d/%s/",year,period.Data());
90   TGridResult *gr = gGrid->Query(path,fileName);
91   Int_t nFiles = gr->GetEntries();
92   printf("================>%d files found\n", nFiles);
93   if (nFiles < 1) return;
94
95
96   if (nFiles > 1){
97     for (Int_t iFil = 0; iFil <nFiles ; iFil++) { 
98       TString fileNameLong=Form("%s",gr->GetKey(iFil,"turl"));
99       if(!fileNameLong.Contains(recoPass.Data())) continue;
100       if(!fileNameLong.Contains(qaTrain.Data())) continue;
101       TString runNumber=fileNameLong;
102       runNumber.ReplaceAll(Form("alien:///alice/data/%d/%s/",year,period.Data()),"");
103       runNumber.Remove(9,runNumber.Sizeof());
104    
105       Int_t iRun=atoi(runNumber.Data());
106       if(readRun->TestBitNumber(iRun)){ 
107         printf("Run %d aleady in local ntuple -> skipping it\n",iRun);
108         continue;
109       }
110       if(iRun<firstRun) continue;
111       if(iRun>lastRun) continue;
112
113       TString isMerged=fileNameLong;
114       isMerged.Remove(isMerged.Sizeof()-16); 
115       isMerged.Remove(0,isMerged.Sizeof()-5);
116       if(!isMerged.Contains("QA")) continue;
117       printf("Open File %s  Run %d\n",fileNameLong.Data(),iRun);
118
119
120
121       TFile* f=TFile::Open(fileNameLong.Data());  
122
123       TDirectoryFile* df=(TDirectoryFile*)f->Get("SDD_Performance");
124       if(!df){
125         printf("Run %d SDD_Performance MISSING -> Exit\n",iRun);
126         continue;
127       }
128       TList* l=(TList*)df->Get("coutputRP");
129       if(!df){
130         printf("Run %d coutputRP TList MISSING -> Exit\n",iRun);
131         continue;
132       }  
133       TH1F* hcllay=(TH1F*)l->FindObject("hCluInLay");
134       Float_t fracT[6]={0.,0.,0.,0.,0.,0.};
135       Float_t efracT[6]={0.,0.,0.,0.,0.,0.};
136       if(hcllay->GetBinContent(1)>0){
137         for(Int_t iLay=0; iLay<6; iLay++){
138           fracT[iLay]=hcllay->GetBinContent(iLay+2)/hcllay->GetBinContent(1);
139           efracT[iLay]=TMath::Sqrt(fracT[iLay]*(1-fracT[iLay])/hcllay->GetBinContent(1));
140         }
141       }
142       TH1F* hmodT=(TH1F*)l->FindObject("hTPMod");
143       TH1F* hgamod=(TH1F*)l->FindObject("hGAMod");
144       Int_t bestMod=0;
145       for(Int_t iMod=0; iMod<260;iMod++){
146         Int_t gda=(Int_t)hgamod->GetBinContent(iMod+1);
147         if(gda>bestMod) bestMod=gda;
148       }
149       Int_t nChunks=1;
150       if(bestMod>512){
151         nChunks=(Int_t)(bestMod/512.+0.5);
152       }
153       hgamod->Scale(1./nChunks);
154
155       TH1F* hev=(TH1F*)l->FindObject("hNEvents");
156       Int_t nTotEvents=hev->GetBinContent(2);
157       Int_t nTrigEvents=hev->GetBinContent(3);
158       Int_t nEvents=nTotEvents;
159       printf("Run %d Number of Events = %d Triggered=%d\n",iRun,nTotEvents,nTrigEvents);
160       if(nTrigEvents>0){ 
161         nEvents=nTrigEvents;
162       }
163       if(nTotEvents==0) continue;
164       Int_t nModGood3=0;
165       Int_t nModGood4=0;
166       Int_t nModBadAn=0;
167       Float_t sumtp3=0;
168       Float_t sumtp4=0;
169       Float_t sumEtp3=0;
170       Float_t sumEtp4=0;
171       for(Int_t iMod=0; iMod<260; iMod++){
172         Float_t tps=hmodT->GetBinContent(iMod+1);
173         Float_t ga=hgamod->GetBinContent(iMod+1);
174         if(ga<500) nModBadAn++;
175         Float_t tpsN=0.;
176         Float_t etpsN=0.;
177         if(ga>0){
178           tpsN=tps/ga/(Float_t)nEvents;
179           etpsN=TMath::Sqrt(tps)/ga/(Float_t)nEvents;
180           if(iMod<84){
181             sumtp3+=tpsN;
182             sumEtp3+=(etpsN*etpsN);
183             nModGood3++;
184           }
185           else{
186             sumtp4+=tpsN;
187             sumEtp4+=(etpsN*etpsN);
188             nModGood4++;
189           }
190         }
191       }
192
193       TH1F* hapmod=(TH1F*)l->FindObject("hAllPmod");
194       TH1F* hgpmod=(TH1F*)l->FindObject("hGoodPmod");
195       //     TH1F* hmpmod=(TH1F*)l->FindObject("hMissPmod");
196       TH1F* hbrmod=(TH1F*)l->FindObject("hBadRegmod");
197       TH1F* hskmod=(TH1F*)l->FindObject("hSkippedmod");
198       TH1F* hoamod=(TH1F*)l->FindObject("hOutAccmod");
199       TH1F* hnrmod=(TH1F*)l->FindObject("hNoRefitmod");
200       Int_t nBelow95=0;
201       Int_t nBelow80=0;
202       Int_t nBelow60=0;
203       Int_t nZeroP=0;
204       for(Int_t imod=0; imod<260;imod++){
205         Float_t numer=hgpmod->GetBinContent(imod+1)+hbrmod->GetBinContent(imod+1)+hoamod->GetBinContent(imod+1)+hnrmod->GetBinContent(imod+1)+hskmod->GetBinContent(imod+1);
206         Float_t denom=hapmod->GetBinContent(imod+1);
207         if(denom>0){
208           Float_t eff=numer/denom;
209           if(eff<0.95) nBelow95++;
210           if(eff<0.80) nBelow80++;
211           if(eff<0.60) nBelow60++;
212         }
213         if(hmodT->GetBinContent(imod+1)<1.){
214           nZeroP++;
215         }       
216       }
217
218       TH1F* htimT=(TH1F*)l->FindObject("hDrTimTPAll");
219       TH1F* htimTe=(TH1F*)l->FindObject("hDrTimTPExtra");
220       
221       Double_t fracExtra=0.;
222       Double_t errFracExtra=0.;
223       if(htimT->GetEntries()>0){
224         fracExtra=htimTe->GetEntries()/htimT->GetEntries();
225         errFracExtra=TMath::Sqrt(htimTe->GetEntries())/htimT->GetEntries();
226       }
227       Double_t averPoints=0.;
228       Double_t cntBins=0.;
229       for(Int_t iBin=1; iBin<=htimT->GetNbinsX(); iBin++){
230         Float_t tim=htimT->GetBinCenter(iBin);
231         if(tim>2000. && tim<4000.){
232           averPoints+=htimT->GetBinContent(iBin);
233           cntBins+=1;
234         }
235       }
236       Double_t minTime=-999.;
237       Double_t errMinTime=0.;
238       if(cntBins>0){ 
239         averPoints/=cntBins;
240         for(Int_t iBin=1; iBin<=htimT->GetNbinsX(); iBin++){
241           if(htimT->GetBinContent(iBin)>0.5*averPoints){
242             minTime=htimT->GetBinCenter(iBin);
243             errMinTime=0.5*htimT->GetBinWidth(iBin);
244             break;
245           }
246         }
247       }
248       TH1F* hSigTim0=(TH1F*)l->FindObject("hSigTimeInt0");
249       TH1F* hSigTim5=(TH1F*)l->FindObject("hSigTimeInt5");
250
251       Int_t index=0;
252       xnt[index++]=(Float_t)iRun;
253       xnt[index++]=fracT[0];
254       xnt[index++]=efracT[0];
255       xnt[index++]=fracT[1];
256       xnt[index++]=efracT[1];
257       xnt[index++]=fracT[2];
258       xnt[index++]=efracT[2];
259       xnt[index++]=fracT[3];
260       xnt[index++]=efracT[3];
261       xnt[index++]=fracT[4];
262       xnt[index++]=efracT[4];
263       xnt[index++]=fracT[5];
264       xnt[index++]=efracT[5];
265       xnt[index++]=sumtp3/nModGood3;
266       xnt[index++]=TMath::Sqrt(sumEtp3)/nModGood3;
267       xnt[index++]=sumtp4/nModGood4;
268       xnt[index++]=TMath::Sqrt(sumEtp4)/nModGood4;
269       xnt[index++]=minTime;
270       xnt[index++]=errMinTime;
271       xnt[index++]=htimT->GetMean();
272       xnt[index++]=htimT->GetMeanError();
273       xnt[index++]=fracExtra;
274       xnt[index++]=errFracExtra;
275       xnt[index++]=hSigTim0->GetMean();
276       xnt[index++]=hSigTim0->GetMeanError();
277       xnt[index++]=hSigTim5->GetMean();
278       xnt[index++]=hSigTim5->GetMeanError();
279       xnt[index++]=(Float_t)nBelow95;
280       xnt[index++]=(Float_t)nBelow80;
281       xnt[index++]=(Float_t)nBelow60;
282       xnt[index++]=(Float_t)nZeroP;
283       ntsdd->Fill(xnt);
284     }
285   }
286
287   TFile* outfil=new TFile(outFilNam.Data(),"recreate");
288   outfil->cd();
289   ntsdd->Write();
290   outfil->Close();
291
292   MakePlots(outFilNam);
293
294 }
295
296 void MakePlots(TString ntupleFileName){
297   TFile* fil=new TFile(ntupleFileName.Data(),"read");
298   if(!fil){
299     printf("File with ntuple does not exist\n");
300     return;
301   }
302   TNtuple* ntsdd=(TNtuple*)fil->Get("ntsdd");
303
304   Float_t nrun;
305   Float_t meanTrPts3,errmeanTrPts3,meanTrPts4,errmeanTrPts4;
306   Float_t minDrTime,errminDrTime,meanDrTime,errmeanDrTime;
307   Float_t fracTrackWithClu1,fracTrackWithClu2,errfracTrackWithClu1,errfracTrackWithClu2;
308   Float_t fracTrackWithClu3,fracTrackWithClu4,errfracTrackWithClu3,errfracTrackWithClu4;
309   Float_t fracTrackWithClu5,fracTrackWithClu6,errfracTrackWithClu5,errfracTrackWithClu6;
310   Float_t fracExtra,errfracExtra;
311   Float_t meandEdxTB0,errmeandEdxTB0,meandEdxTB5,errmeandEdxTB5;
312   Float_t nMod95,nMod80,nMod60,nModEmpty;
313   
314   ntsdd->SetBranchAddress("nrun",&nrun);
315   ntsdd->SetBranchAddress("fracTrackWithClu1",&fracTrackWithClu1);
316   ntsdd->SetBranchAddress("errfracTrackWithClu1",&errfracTrackWithClu1);
317   ntsdd->SetBranchAddress("fracTrackWithClu2",&fracTrackWithClu2);
318   ntsdd->SetBranchAddress("errfracTrackWithClu2",&errfracTrackWithClu2);
319   ntsdd->SetBranchAddress("fracTrackWithClu3",&fracTrackWithClu3);
320   ntsdd->SetBranchAddress("errfracTrackWithClu3",&errfracTrackWithClu3);
321   ntsdd->SetBranchAddress("fracTrackWithClu4",&fracTrackWithClu4);
322   ntsdd->SetBranchAddress("errfracTrackWithClu4",&errfracTrackWithClu4);
323   ntsdd->SetBranchAddress("fracTrackWithClu5",&fracTrackWithClu5);
324   ntsdd->SetBranchAddress("errfracTrackWithClu5",&errfracTrackWithClu5);
325   ntsdd->SetBranchAddress("fracTrackWithClu6",&fracTrackWithClu6);
326   ntsdd->SetBranchAddress("errfracTrackWithClu6",&errfracTrackWithClu6);
327   ntsdd->SetBranchAddress("nMod95",&nMod95);
328   ntsdd->SetBranchAddress("nMod80",&nMod80);
329   ntsdd->SetBranchAddress("nMod60",&nMod60);
330   ntsdd->SetBranchAddress("nModEmpty",&nModEmpty);
331
332   ntsdd->SetBranchAddress("meanTrPts3",&meanTrPts3);
333   ntsdd->SetBranchAddress("errmeanTrPts3",&errmeanTrPts3);
334   ntsdd->SetBranchAddress("meanTrPts4",&meanTrPts4);
335   ntsdd->SetBranchAddress("errmeanTrPts4",&errmeanTrPts4);
336   ntsdd->SetBranchAddress("minDrTime",&minDrTime);
337   ntsdd->SetBranchAddress("errminDrTime",&errminDrTime);
338   ntsdd->SetBranchAddress("meanDrTime",&meanDrTime);
339   ntsdd->SetBranchAddress("errmeanDrTime",&errmeanDrTime);
340   ntsdd->SetBranchAddress("fracExtra",&fracExtra);
341   ntsdd->SetBranchAddress("errfracExtra",&errfracExtra);
342   ntsdd->SetBranchAddress("meandEdxTB0",&meandEdxTB0);
343   ntsdd->SetBranchAddress("errmeandEdxTB0",&errmeandEdxTB0);
344   ntsdd->SetBranchAddress("meandEdxTB5",&meandEdxTB5);
345   ntsdd->SetBranchAddress("errmeandEdxTB5",&errmeandEdxTB5);
346
347   TH1F* histotrp3=new TH1F("histotrp3","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
348   TH1F* histotrp4=new TH1F("histotrp4","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
349   TH1F* histominTime=new TH1F("histominTime","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
350   TH1F* histomeanTime=new TH1F("histomeanTime","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
351   TH1F* histofracExtra=new TH1F("histofracExtra","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
352   TH1F* histodEdxTB0=new TH1F("histodEdxTB0","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
353   TH1F* histodEdxTB5=new TH1F("histodEdxTB5","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
354   TH1F* histoTrackClu1=new TH1F("histoTrackClu1","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
355   TH1F* histoTrackClu2=new TH1F("histoTrackClu2","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
356   TH1F* histoTrackClu3=new TH1F("histoTrackClu3","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
357   TH1F* histoTrackClu4=new TH1F("histoTrackClu4","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
358   TH1F* histoTrackClu5=new TH1F("histoTrackClu5","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
359   TH1F* histoTrackClu6=new TH1F("histoTrackClu6","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
360
361   TH1F* histoNmodEffBelow95=new TH1F("histoNmodEffBelow95","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
362   TH1F* histoNmodEffBelow80=new TH1F("histoNmodEffBelow80","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
363   TH1F* histoNmodEffBelow60=new TH1F("histoNmodEffBelow60","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
364   TH1F* histoNmodEmpty=new TH1F("histoNmodEmpty","",(Int_t)ntsdd->GetEntries(),0.,ntsdd->GetEntries());
365
366   for(Int_t i=0; i<ntsdd->GetEntries();i++){
367     ntsdd->GetEvent(i);
368     histoTrackClu1->SetBinContent(i+1,fracTrackWithClu1);
369     histoTrackClu1->SetBinError(i+1,errfracTrackWithClu1);
370     histoTrackClu1->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
371     histoTrackClu2->SetBinContent(i+1,fracTrackWithClu2);
372     histoTrackClu2->SetBinError(i+1,errfracTrackWithClu2);
373     histoTrackClu2->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
374     histoTrackClu3->SetBinContent(i+1,fracTrackWithClu3);
375     histoTrackClu3->SetBinError(i+1,errfracTrackWithClu3);
376     histoTrackClu3->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
377     histoTrackClu4->SetBinContent(i+1,fracTrackWithClu4);
378     histoTrackClu4->SetBinError(i+1,errfracTrackWithClu4);
379     histoTrackClu4->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
380     histoTrackClu5->SetBinContent(i+1,fracTrackWithClu5);
381     histoTrackClu5->SetBinError(i+1,errfracTrackWithClu5);
382     histoTrackClu5->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
383     histoTrackClu6->SetBinContent(i+1,fracTrackWithClu6);
384     histoTrackClu6->SetBinError(i+1,errfracTrackWithClu6);
385     histoTrackClu6->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
386     histominTime->SetBinContent(i+1,minDrTime);
387     histominTime->SetBinError(i+1,errminDrTime);
388     histominTime->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
389     histomeanTime->SetBinContent(i+1,meanDrTime);
390     histomeanTime->SetBinError(i+1,errmeanDrTime);
391     histomeanTime->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
392     histotrp3->SetBinContent(i+1,meanTrPts3);
393     histotrp3->SetBinError(i+1,errmeanTrPts3);
394     histotrp3->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
395     histotrp4->SetBinContent(i+1,meanTrPts4);
396     histotrp4->SetBinError(i+1,errmeanTrPts3);
397     histotrp4->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
398     histofracExtra->SetBinContent(i+1,fracExtra);
399     histofracExtra->SetBinError(i+1,errfracExtra);
400     histofracExtra->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
401     histodEdxTB0->SetBinContent(i+1,meandEdxTB0);
402     histodEdxTB0->SetBinError(i+1,errmeandEdxTB0);
403     histodEdxTB0->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
404     histodEdxTB5->SetBinContent(i+1,meandEdxTB5);
405     histodEdxTB5->SetBinError(i+1,errmeandEdxTB5);
406     histodEdxTB5->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
407
408     histoNmodEffBelow95->SetBinContent(i+1,nMod95);
409     histoNmodEffBelow95->SetBinError(i+1,0.0001);
410     histoNmodEffBelow95->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
411     histoNmodEffBelow80->SetBinContent(i+1,nMod80);
412     histoNmodEffBelow80->SetBinError(i+1,0.0001);
413     histoNmodEffBelow80->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
414     histoNmodEffBelow60->SetBinContent(i+1,nMod60);
415     histoNmodEffBelow60->SetBinError(i+1,0.0001);
416     histoNmodEffBelow60->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
417     histoNmodEmpty->SetBinContent(i+1,nModEmpty);
418     histoNmodEmpty->SetBinError(i+1,0.0001);
419     histoNmodEmpty->GetXaxis()->SetBinLabel(i+1,Form("%d",(Int_t)nrun));
420   }
421
422   gStyle->SetOptStat(0);
423
424
425   TCanvas* c1=new TCanvas("c1","Occupancy");
426   histotrp3->SetLineColor(1);
427   histotrp3->SetMarkerStyle(20);
428   histotrp3->Draw();
429   histotrp3->SetMinimum(0.);
430   histotrp4->SetLineColor(2);
431   histotrp4->SetMarkerColor(2);
432   histotrp4->SetMarkerStyle(22);
433   histotrp4->Draw("same");
434   histotrp3->GetYaxis()->SetTitle("Track Point Occupancy");
435   TLegend* leg=new TLegend(0.7,0.15,0.88,0.35);
436   TLegendEntry* ent=leg->AddEntry(histotrp3,"Layer3","PL");
437   ent=leg->AddEntry(histotrp4,"Layer4","PL");
438   ent->SetTextColor(histotrp4->GetMarkerColor());
439   leg->SetFillStyle(0);
440   leg->Draw();
441   c1->Update();
442
443   TCanvas* c2=new TCanvas("c2","DriftTime",600,900);
444   c2->Divide(1,2);
445   c2->cd(1);
446   histominTime->Draw();
447   histominTime->SetMinimum(250);
448   histominTime->SetMaximum(1000);
449   histominTime->GetYaxis()->SetTitle("Minimum Drift Time (ns)");
450   c2->cd(2);
451   histomeanTime->Draw();
452   histomeanTime->GetYaxis()->SetTitle("Average Drift Time (ns)");
453   c2->Update();
454
455   TCanvas* c3=new TCanvas("c3","ExtraClusters");
456   histofracExtra->Draw();
457   histofracExtra->GetYaxis()->SetTitle("Fraction of Extra Clusters");
458   c3->Update();
459
460
461   TCanvas* c4=new TCanvas("c4","Charge");
462   histodEdxTB0->SetLineColor(1);
463   histodEdxTB0->SetMarkerStyle(20);
464   histodEdxTB0->Draw();
465   histodEdxTB0->SetMinimum(0.);
466   histodEdxTB5->SetLineColor(4);
467   histodEdxTB5->SetMarkerColor(4);
468   histodEdxTB5->SetMarkerStyle(23);
469   histodEdxTB5->Draw("same");
470   histodEdxTB0->GetYaxis()->SetTitle("<dE/dx> (keV/300 #mum)");
471   TLegend* leg2=new TLegend(0.6,0.15,0.88,0.35);
472   ent=leg2->AddEntry(histodEdxTB0,"Small drift time","PL");
473   ent=leg2->AddEntry(histodEdxTB5,"Large drift time","PL");
474   ent->SetTextColor(histodEdxTB5->GetMarkerColor());
475   leg2->SetFillStyle(0);
476   leg2->Draw();
477   c4->Update();
478
479   TCanvas* c5=new TCanvas("c5","TrackWithSDD");
480   histoTrackClu3->Draw();
481   histoTrackClu3->SetLineColor(1);
482   histoTrackClu3->SetMarkerStyle(20);
483   histoTrackClu3->Draw();
484   histoTrackClu3->SetMinimum(0.);
485   histoTrackClu4->SetLineColor(2);
486   histoTrackClu4->SetMarkerColor(2);
487   histoTrackClu4->SetMarkerStyle(22);
488   histoTrackClu4->Draw("same");
489   histoTrackClu1->SetLineColor(kGray+1);
490   histoTrackClu1->SetMarkerColor(kGray+1);
491   histoTrackClu1->SetMarkerStyle(24);
492   histoTrackClu1->Draw("same");
493   histoTrackClu2->SetLineColor(kGray+2);
494   histoTrackClu2->SetMarkerColor(kGray+2);
495   histoTrackClu2->SetMarkerStyle(26);
496   histoTrackClu2->Draw("same");
497   histoTrackClu5->SetLineColor(4);
498   histoTrackClu5->SetMarkerColor(4);
499   histoTrackClu5->SetMarkerStyle(29);
500   histoTrackClu5->Draw("same");
501   histoTrackClu6->SetLineColor(kBlue+1);
502   histoTrackClu6->SetMarkerColor(kBlue+1);
503   histoTrackClu6->SetMarkerStyle(30);
504   histoTrackClu6->Draw("same");
505   histoTrackClu3->GetYaxis()->SetTitle("Fraction of Tracks with Point in SDD Layers");
506   TLegend* leg3=new TLegend(0.7,0.15,0.88,0.35);
507   ent=leg3->AddEntry(histoTrackClu1,"Layer1","PL");
508   ent->SetTextColor(histoTrackClu1->GetMarkerColor());
509   ent=leg3->AddEntry(histoTrackClu2,"Layer2","PL");
510   ent->SetTextColor(histoTrackClu2->GetMarkerColor());
511   ent=leg3->AddEntry(histoTrackClu3,"Layer3","PL");
512   ent->SetTextColor(histoTrackClu3->GetMarkerColor());
513   ent=leg3->AddEntry(histoTrackClu4,"Layer4","PL");
514   ent->SetTextColor(histoTrackClu4->GetMarkerColor());
515   ent=leg3->AddEntry(histoTrackClu5,"Layer5","PL");
516   ent->SetTextColor(histoTrackClu5->GetMarkerColor());
517   ent=leg3->AddEntry(histoTrackClu6,"Layer6","PL");
518   ent->SetTextColor(histoTrackClu6->GetMarkerColor());
519
520   leg3->SetFillStyle(0);
521   leg3->Draw();
522   c5->Update();
523
524   TCanvas* c6=new TCanvas("c6","Modules with low eff",800,1000);
525   c6->Divide(1,4);
526   c6->cd(1);
527   histoNmodEffBelow95->Draw("E");
528   c6->cd(2);
529   histoNmodEffBelow80->Draw("E");
530   c6->cd(3);
531   histoNmodEffBelow60->Draw("E");
532   c6->cd(4);
533   histoNmodEmpty->Draw("E");
534
535 }
536
537 Double_t LangausFun(Double_t *x, Double_t *par) {
538
539   //Fit parameters:
540   //par[0]=Width (scale) parameter of Landau density
541   //par[1]=Most Probable (MP, location) parameter of Landau density
542   //par[2]=Total area (integral -inf to inf, normalization constant)
543   //par[3]=Width (sigma) of convoluted Gaussian function
544   //
545   //In the Landau distribution (represented by the CERNLIB approximation), 
546   //the maximum is located at x=-0.22278298 with the location parameter=0.
547   //This shift is corrected within this function, so that the actual
548   //maximum is identical to the MP parameter.
549
550   // Numeric constants
551   Double_t invsq2pi = 0.3989422804014;   // (2 pi)^(-1/2)
552   Double_t mpshift  = -0.22278298;       // Landau maximum location
553
554   // Control constants
555   Double_t np = 100.0;      // number of convolution steps
556   Double_t sc =   5.0;      // convolution extends to +-sc Gaussian sigmas
557
558   // Variables
559   Double_t xx;
560   Double_t mpc;
561   Double_t fland;
562   Double_t sum = 0.0;
563   Double_t xlow,xupp;
564   Double_t step;
565   Double_t i;
566
567
568   // MP shift correction
569   mpc = par[1] - mpshift * par[0]; 
570
571   // Range of convolution integral
572   xlow = x[0] - sc * par[3];
573   xupp = x[0] + sc * par[3];
574
575   step = (xupp-xlow) / np;
576
577   // Convolution integral of Landau and Gaussian by sum
578   for(i=1.0; i<=np/2; i++) {
579     xx = xlow + (i-.5) * step;
580     fland = TMath::Landau(xx,mpc,par[0]) / par[0];
581     sum += fland * TMath::Gaus(x[0],xx,par[3]);
582
583     xx = xupp - (i-.5) * step;
584     fland = TMath::Landau(xx,mpc,par[0]) / par[0];
585     sum += fland * TMath::Gaus(x[0],xx,par[3]);
586   }
587
588   return (par[2] * step * sum * invsq2pi / par[3]);
589 }