]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/macrosSDD/ShowResponseSDD.C
Minor updates in two SDD macros
[u/mrichter/AliRoot.git] / ITS / macrosSDD / ShowResponseSDD.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TFile.h>
3 #include <TH1F.h>
4 #include <TH2I.h>
5 #include <TGraph.h>
6 #include <TExec.h>
7 #include <TStyle.h>
8 #include <TString.h>
9 #include <TSystem.h>
10 #include <TGrid.h>
11 #include <TLatex.h>
12 #include <TCanvas.h>
13 #include <TObjArray.h>
14 #include <TPaveStats.h>
15 #include "AliCDBEntry.h"
16 #include "AliITSresponseSDD.h"
17 #endif
18
19 /*  $Id: ShowResponseSDD.C 44072 2010-10-04 15:48:32Z prino $    */
20
21 // Macro to plot the calibration parameters 
22 // (time zero, ADCtokeV and Vdrift correction)
23 // from the OCDB file created from SDD offline calibration 
24 // (OCDB/ITS/Calib/RespSDD)
25
26 void ShowResponseSDD(TString filename="$ALICE_ROOT/OCDB/ITS/Calib/RespSDD/Run0_999999999_v0_s0.root"){
27   if(filename.Contains("alien")){
28     TGrid::Connect("alien:");
29   }
30   TFile* fr=TFile::Open(filename.Data());
31   AliCDBEntry* e=(AliCDBEntry*)fr->Get("AliCDBEntry");
32   e->PrintMetaData();
33   AliITSresponseSDD* r=(AliITSresponseSDD*)e->GetObject();
34   TH1F* hTimeZero=new TH1F("hTimeZero","",260,239.5,499.5);
35   TH1F* hVdriftCorrLeft=new TH1F("hVdriftCorrLeft","",260,239.5,499.5);
36   TH1F* hVdriftCorrRight=new TH1F("hVdriftCorrRight","",260,239.5,499.5);
37   TH1F* hdistVdriftCorrLeft=new TH1F("hdistVdriftCorrLeft","",100,-0.2,0.2);
38   TH1F* hdistVdriftCorrRight=new TH1F("hdistVdriftCorrRight","",100,-0.2,0.2);
39   TH1F* hADCtokeV=new TH1F("hADCtokeV","",260,239.5,499.5);
40   TH1F* hADCvsTime=new TH1F("hADCvsTime","",260,239.5,499.5);
41   Float_t averTz=0.;
42   Float_t averCv0=0.;
43   Float_t averCv1=0.;
44   Float_t averAk=0.;
45   Float_t averAt=0.;
46   for(Int_t iMod=240; iMod<500; iMod++){
47     Float_t tz=r->GetTimeZero(iMod);
48     Float_t cv0=r->GetDeltaVDrift(iMod,kFALSE);
49     Float_t cv1=r->GetDeltaVDrift(iMod,kTRUE);
50     Float_t ak=r->GetADCtokeV(iMod);
51     Float_t at=r->GetADCvsDriftTime(iMod);
52     hTimeZero->SetBinContent(iMod-240+1,tz);    
53     hVdriftCorrLeft->SetBinContent(iMod-240+1,cv0);
54     hVdriftCorrRight->SetBinContent(iMod-240+1,cv1);
55     hdistVdriftCorrLeft->Fill(cv0);
56     hdistVdriftCorrRight->Fill(cv1);
57     hADCtokeV->SetBinContent(iMod-240+1,ak);
58     hADCvsTime->SetBinContent(iMod-240+1,at);
59     averTz+=tz;
60     averCv0+=cv0;
61     averCv1+=cv1;
62     averAk+=ak;
63     averAt+=at;
64   }
65   averTz/=260.;
66   averCv0/=260.;
67   averCv1/=260.;
68   averAk/=260.;
69   averAt/=260.;
70
71   hTimeZero->SetMarkerStyle(20);
72   hADCtokeV->SetMarkerStyle(20);
73   hADCvsTime->SetMarkerStyle(20);
74   hVdriftCorrLeft->SetMarkerStyle(22);
75   hVdriftCorrLeft->SetMarkerColor(2);
76   hVdriftCorrRight->SetMarkerStyle(24);
77   hVdriftCorrRight->SetMarkerColor(4);
78   hTimeZero->SetMaximum(hTimeZero->GetMaximum()*1.2);
79   hADCtokeV->SetMaximum(hADCtokeV->GetMaximum()*1.2);
80   hADCvsTime->SetMaximum(hADCvsTime->GetMaximum()*1.2);
81   Float_t maxV=TMath::Max(hVdriftCorrLeft->GetMaximum(),hVdriftCorrRight->GetMaximum());
82   Float_t minV=TMath::Min(hVdriftCorrLeft->GetMinimum(),hVdriftCorrRight->GetMinimum());
83   Float_t scale=TMath::Max(TMath::Abs(maxV),TMath::Abs(minV));
84   if(scale<0.02){
85     hVdriftCorrLeft->SetMinimum(-0.05);
86     hVdriftCorrLeft->SetMaximum(0.05);
87   }else{
88     hVdriftCorrLeft->SetMaximum(1.2*scale);
89     hVdriftCorrLeft->SetMinimum(-1.2*scale);
90   }
91
92   //  gStyle->SetOptStat(0);
93
94   printf("Charge vs. Time correction factor = %f\n",r->GetChargevsTime());
95
96   TCanvas* c1=new TCanvas("c1","Time Zero");
97   hTimeZero->SetStats(0);
98   hTimeZero->Draw("P");
99   hTimeZero->GetXaxis()->SetTitle("Module Id");
100   hTimeZero->GetYaxis()->SetTitle("Time Zero (ns)");
101   TLatex* tt=new TLatex(0.2,0.8,Form("Average Tzero = %.2f",averTz));
102   tt->SetNDC();
103   tt->Draw();
104   c1->Modified();
105
106   TCanvas* c2=new TCanvas("c2","Vdrift Corr",800,900);
107   c2->Divide(1,2);
108   c2->cd(1);
109   hVdriftCorrLeft->SetStats(0);
110   hVdriftCorrLeft->Draw("P");
111   hVdriftCorrRight->Draw("PSAME");
112   hVdriftCorrLeft->GetXaxis()->SetTitle("Module Id");
113   hVdriftCorrLeft->GetYaxis()->SetTitle("Vdrift correction (#mum/ns)");
114   TLatex* tc0=new TLatex(0.2,0.8,Form("Average Vdrift corr Left = %.3f",averCv0));
115   tc0->SetNDC();
116   tc0->SetTextColor(2);
117   tc0->Draw();
118   TLatex* tc1=new TLatex(0.2,0.72,Form("Average Vdrift corr Right = %.3f",averCv1));
119   tc1->SetNDC();
120   tc1->SetTextColor(4);
121   tc1->Draw();
122   c2->cd(2);
123   hdistVdriftCorrLeft->SetLineColor(2);
124   hdistVdriftCorrLeft->GetXaxis()->SetTitle("Vdrift correction (#mum/ns)");
125   hdistVdriftCorrLeft->Draw(); 
126   c2->Update();
127   TPaveStats *stL=(TPaveStats*)hdistVdriftCorrLeft->GetListOfFunctions()->FindObject("stats");
128   stL->SetY1NDC(0.71);
129   stL->SetY2NDC(0.9);
130   stL->SetTextColor(2);
131   hdistVdriftCorrRight->SetLineColor(4);
132   hdistVdriftCorrRight->Draw("SAMES");
133   c2->Update();
134   TPaveStats *stR=(TPaveStats*)hdistVdriftCorrRight->GetListOfFunctions()->FindObject("stats");
135   stR->SetY1NDC(0.51);
136   stR->SetY2NDC(0.7);
137   stR->SetTextColor(4);
138
139   TCanvas* c3=new TCanvas("c3","ADC calib",800,900);
140   c3->Divide(1,2);
141   c3->cd(1);
142   hADCvsTime->SetStats(0);
143   hADCvsTime->Draw("P");
144   hADCvsTime->GetXaxis()->SetTitle("Module Id");
145   hADCvsTime->GetYaxis()->SetTitle("ADC vs. Drift Time Slope (ADC/ns)");
146   TLatex* ts=new TLatex(0.2,0.8,Form("Average ADCvsTime = %.3f",averAt));
147   ts->SetNDC();
148   ts->Draw(); 
149   c3->cd(2);
150   hADCtokeV->SetStats(0);
151   hADCtokeV->Draw("P");
152   hADCtokeV->GetXaxis()->SetTitle("Module Id");
153   hADCtokeV->GetYaxis()->SetTitle("ADC to keV conv. factor");
154   TLatex* ta=new TLatex(0.2,0.8,Form("Average ADCtokeV = %.3f",averAk));
155   ta->SetNDC();
156   ta->Draw(); 
157
158 }