]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/macrosSDD/ShowCorrMapsSDD.C
New plot for trending dE/dx on Lay3 and Lay4 separately
[u/mrichter/AliRoot.git] / ITS / macrosSDD / ShowCorrMapsSDD.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TFile.h>
3 #include <TH1F.h>
4 #include <TH2F.h>
5 #include <TStyle.h>
6 #include <TSystem.h>
7 #include <TString.h>
8 #include <TGrid.h>
9 #include <TCanvas.h>
10 #include <TObjArray.h>
11 #include "AliCDBEntry.h"
12 #include "AliITSCorrMapSDD.h"
13 #include "AliITSMapSDD.h"
14 #endif
15
16
17 // Macro to plot the SDD correction maps 
18 // from the OCDB file (OCDB/ITS/Calib/MapsTimeSDD)
19 //
20 // Origin: F. Prino (prino@to.infn.it)
21
22 void ShowCorrMapsSDD(TString filname="alien:///alice/data/2010/OCDB/ITS/Calib/MapsTimeSDD/Run117224_999999999_v3_s0.root"){
23
24   if(filname.Contains("alien")){
25     TGrid::Connect("alien:");
26   }
27   TFile *fil =TFile::Open(filname.Data());
28     
29   AliCDBEntry* ent=(AliCDBEntry*)fil->Get("AliCDBEntry");
30   TObjArray *maptSDD = (TObjArray *)ent->GetObject();
31   printf("Entries in array=%d\n",maptSDD->GetEntriesFast());
32   TString psnm0 = "mapsSDD.ps[";
33   TString psnm1 = "mapsSDD.ps";
34   TString psnm2 = "mapsSDD.ps]";
35   Bool_t oldMapFormat=kFALSE;
36   TObject* objmap=(TObject*)maptSDD->At(0);
37   TString cname(objmap->ClassName());
38   if(cname.CompareTo("AliITSMapSDD")==0){ 
39     oldMapFormat=kTRUE;
40     printf("SDD Maps converted to new format\n");
41   }
42  
43
44   TCanvas * c3=new TCanvas("c3","Layer 3",1200,900);
45   c3->Print(psnm0.Data());
46   gStyle->SetPalette(1);
47   gStyle->SetOptStat(0);
48   Int_t cntpad=0;
49   for(Int_t i=0; i<84; i++){
50     if(i%6==0){
51       c3->cd();
52       c3->Modified();
53       c3->Update();
54       if(i) c3->Print(psnm1.Data());
55       c3->Clear();
56       c3->Divide(4,3);
57       cntpad=0;
58     }
59     Int_t index0=i*2;
60     Int_t index1=i*2+1;
61     AliITSCorrMapSDD* map0 = 0;
62     AliITSCorrMapSDD* map1 = 0;
63     if(oldMapFormat){ 
64       AliITSMapSDD* oldmap0=(AliITSMapSDD*)maptSDD->At(index0);
65       AliITSMapSDD* oldmap1=(AliITSMapSDD*)maptSDD->At(index1);
66       map0=oldmap0->ConvertToNewFormat();
67       map1=oldmap1->ConvertToNewFormat();
68     }else{
69       map0=(AliITSCorrMapSDD*)maptSDD->At(index0);
70       map1=(AliITSCorrMapSDD*)maptSDD->At(index1);
71     }
72
73     printf("Module %s Entries in map %dx%d\n",map0->GetName(),map0->GetNBinsAnode(),map0->GetNBinsDrift());
74     if(map0->GetNBinsAnode()==1){
75       TH1F* hp0=map0->GetMapProfile();
76       TH1F* hp1=map1->GetMapProfile();
77       if(hp0->GetMinimum()>-0.000001) hp0->SetMinimum(-10);
78       if(hp0->GetMaximum()<0.000001) hp0->SetMaximum(10);
79       if(hp1->GetMinimum()>-0.000001) hp1->SetMinimum(-10);
80       if(hp1->GetMaximum()<0.000001) hp1->SetMaximum(10);
81       hp0->SetTitle(Form("Module %d - Left",i+240));
82       hp1->SetTitle(Form("Module %d - Right",i+240));
83       hp0->GetXaxis()->SetTitle("Drift distance (mm)");
84       hp0->GetYaxis()->SetTitle("Correction (#mum)");
85       hp1->GetXaxis()->SetTitle("Drift distance (mm)");
86       hp1->GetYaxis()->SetTitle("Correction (#mum)");
87       c3->cd(++cntpad);
88       hp0->Draw();
89       c3->cd(++cntpad);
90       hp1->Draw();
91     }else{
92       TH2F* hmap0=map0->GetMapHisto();
93       TH2F* hmap1=map1->GetMapHisto();
94       c3->cd(++cntpad);
95       hmap0->Draw("colz");
96       c3->cd(++cntpad);
97       hmap1->Draw("colz");
98     }
99   }
100   c3->cd();
101   c3->Modified();
102   c3->Update();
103   c3->Print(psnm1.Data());
104
105   TCanvas * c4=new TCanvas("c4","Layer 4",1200,900);
106   c4->Divide(4,4);
107   gStyle->SetPalette(1);
108   for(Int_t i=84; i<260; i++){
109     Int_t j=i-84;
110     if(j%8==0){
111       c4->cd();
112       c4->Modified();
113       c4->Update();
114       if(j) c4->Print(psnm1.Data());
115       c4->Clear();
116       c4->Divide(4,4);
117       cntpad=0;
118     }
119     Int_t index0=i*2;
120     Int_t index1=i*2+1;
121     AliITSCorrMapSDD* map0 = 0;
122     AliITSCorrMapSDD* map1 = 0;
123     if(oldMapFormat){ 
124       AliITSMapSDD* oldmap0=(AliITSMapSDD*)maptSDD->At(index0);
125       AliITSMapSDD* oldmap1=(AliITSMapSDD*)maptSDD->At(index1);
126       map0=oldmap0->ConvertToNewFormat();
127       map1=oldmap1->ConvertToNewFormat();
128     }else{
129       map0=(AliITSCorrMapSDD*)maptSDD->At(index0);
130       map1=(AliITSCorrMapSDD*)maptSDD->At(index1);
131     }
132     printf("Module %d Entries in map %dx%d\n",i,map0->GetNBinsAnode(),map0->GetNBinsDrift());
133     if(map0->GetNBinsAnode()==1){
134       TH1F* hp0=map0->GetMapProfile();
135       TH1F* hp1=map1->GetMapProfile();
136       if(hp0->GetMinimum()>-0.000001) hp0->SetMinimum(-10);
137       if(hp0->GetMaximum()<0.000001) hp0->SetMaximum(10);
138       if(hp1->GetMinimum()>-0.000001) hp1->SetMinimum(-10);
139       if(hp1->GetMaximum()<0.000001) hp1->SetMaximum(10);
140       hp0->SetTitle(Form("Module %d - Left",i+240));
141       hp1->SetTitle(Form("Module %d - Right",i+240));
142       hp0->GetXaxis()->SetTitle("Drift distance (mm)");
143       hp0->GetYaxis()->SetTitle("Correction (#mum)");
144       hp1->GetXaxis()->SetTitle("Drift distance (mm)");
145       hp1->GetYaxis()->SetTitle("Correction (#mum)");
146       c4->cd(++cntpad);
147       hp0->Draw();
148       c4->cd(++cntpad);
149       hp1->Draw();
150     }else{
151       TH2F* hmap0=map0->GetMapHisto();
152       TH2F* hmap1=map1->GetMapHisto();
153       c4->cd(++cntpad);
154       hmap0->Draw("colz");
155       c4->cd(++cntpad);
156       hmap1->Draw("colz");
157     }
158   }
159   c4->cd();
160   c4->Modified();
161   c4->Update();
162   c4->Print(psnm1.Data());
163   c4->Print(psnm2.Data());
164
165 }