]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/Upgrade/macros/CheckConsistency_DistCorr.C
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / TPC / Upgrade / macros / CheckConsistency_DistCorr.C
1 #include "TString.h"
2 #include "TFile.h"
3 #include "TCanvas.h"
4 #include "TMath.h"
5 #include "TTree.h"
6 #include "TTreeStream.h"
7 #include "TROOT.h"
8 #include "TSystem.h"
9 #include "TStyle.h"
10 #include "AliTPCCorrection.h"
11 #include "AliTPCCorrectionLookupTable.h"
12 #include <AliToyMCEventGenerator.h>
13
14
15
16 void makeComparisonTree(TString filename, TString addToName)
17 {
18   
19   AliTPCCorrectionLookupTable *fTPCCorrection2=0x0;
20
21   Bool_t doScaling=kTRUE;
22   if (filename.Contains(":")) {
23     TObjArray *arr=filename.Tokenize(":");
24     TString s2(arr->At(1)->GetName());
25     if (s2.Contains("-scale")) {
26       doScaling=kFALSE;
27       s2.ReplaceAll("-scale","");
28     }
29     TFile f2(s2);
30     gROOT->cd();
31     fTPCCorrection2=(AliTPCCorrectionLookupTable*)f2.Get("map");
32     f2.Close();
33     filename=arr->At(0)->GetName();
34     delete arr;
35   }
36   TFile fn(filename.Data());
37   gROOT->cd();
38   AliTPCCorrectionLookupTable *fTPCCorrection=(AliTPCCorrectionLookupTable*)fn.Get("map");
39   fn.Close();
40   if (fTPCCorrection2 && doScaling) {
41     Float_t dummy=0;
42     fTPCCorrection2->SetCorrScaleFactor(AliToyMCEventGenerator::GetSCScalingFactor(fTPCCorrection, fTPCCorrection2,dummy));
43     
44   }
45 //   fTPCCorrection->BuildExactInverse();
46
47 //   TFile f("/tmp/corrTest.Root","recreate");
48 //   fTPCCorrection->Write("map");
49 //   f.Close();
50   
51   TString outFile=addToName;
52   outFile.Append(".root");
53   TTreeSRedirector *sred=new TTreeSRedirector(outFile.Data());
54   
55   Float_t dx[3]={0,0,0};
56   
57   for (Float_t iz=-245; iz<=245; iz+=10) {
58     Short_t roc=(iz>=0)?0:18;
59     for (Float_t ir=86; ir<250; ir+=10) {
60       for (Float_t iphi=0; iphi<TMath::TwoPi(); iphi+=10*TMath::DegToRad()){
61         Float_t x=ir*(Float_t)TMath::Cos(iphi);
62         Float_t y=ir*(Float_t)TMath::Sin(iphi);
63         Float_t x3[3]={x,y,iz};
64         (*sred) << "t" <<
65         "r="   << ir <<
66         "phi=" << iphi <<
67         "x=" << x <<
68         "y=" << y <<
69         "z=" << iz;
70         
71         //distortions
72         fTPCCorrection->GetDistortion(x3,roc,dx);
73         Float_t xd   = x+dx[0];
74         Float_t yd   = y+dx[1];
75         Float_t zd   = iz+dx[2];
76         Float_t rd   = TMath::Sqrt(xd*xd+yd*yd);
77         Float_t phid = TMath::ATan2(yd,xd);
78         if (phid<0) phid+=TMath::TwoPi();
79         (*sred) << "t" <<
80         "xd="   << xd <<
81         "yd="   << yd <<
82         "zd="   << zd <<
83         "rd="   << rd <<
84         "phid=" << phid;
85         
86         // correct back distorted point
87         Float_t xd3[3]={xd,yd,zd};
88
89         if (fTPCCorrection2) {
90           fTPCCorrection2->GetCorrection(xd3,roc,dx);
91         } else {
92           fTPCCorrection->GetCorrection(xd3,roc,dx);
93         }
94         Float_t xdc   = xd+dx[0];
95         Float_t ydc   = yd+dx[1];
96         Float_t zdc   = zd+dx[2];
97         Float_t rdc   = TMath::Sqrt(xdc*xdc+ydc*ydc);
98         Float_t phidc = TMath::ATan2(ydc,xdc);
99         if (phidc<0) phidc+=TMath::TwoPi();
100         
101         (*sred)  << "t" <<
102         "xdc="   << xdc <<
103         "ydc="   << ydc <<
104         "zdc="   << zdc <<
105         "rdc="   << rdc <<
106         "phidc=" << phidc;
107         
108         // write current point
109         (*sred) << "t" <<
110         "\n";
111         
112       }
113     }
114   }
115   
116   delete sred;
117 }
118
119 void makeAllComparisonTrees()
120 {
121   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2_eps5_50kHz_precal.lookup.root","LUT_05");
122   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2_eps10_50kHz_precal.lookup.root","LUT_10");
123   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2_eps20_50kHz_precal.lookup.root","LUT_20");
124 }
125
126 void makeAllComparisonTreesNew()
127 {
128   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps5_50kHz_precal.lookup.root","LUT_05");
129   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps10_50kHz_precal.lookup.root","LUT_10");
130   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_20");
131   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_25");
132   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_30");
133   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_35");
134   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_40");
135 }
136
137 void makeAllComparisonTreesOld()
138 {
139   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps5_50kHz_precal.lookup.root","LUT_05");
140   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps10_50kHz_precal.lookup.root","LUT_10");
141   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_20");
142   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_25");
143   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_30");
144   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_35");
145   makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_40");
146 }
147
148 TCanvas *GetCanvas(TString addToName);
149
150 void makeHistos(TString addToName) {
151   TString filename; //("test_");
152   filename.Append(addToName.Data());
153   filename.Append(".root");
154   TFile f(filename.Data());
155   gROOT->cd();
156   TTree *t=(TTree*)f.Get("t");
157   gStyle->SetTitleX(0.18);
158   gStyle->SetTitleW(1-.18-.1);
159
160   t->SetMarkerStyle(20);
161   t->SetMarkerSize(.8);
162
163   TCanvas *c=0x0;
164   c=GetCanvas(addToName+"_zRes");
165   t->Draw("zdc-z:z:r","","colz");
166   c->SaveAs(Form("%s_zRes.png",addToName.Data()));
167   //
168   c=GetCanvas(addToName+"_rRes");
169   t->Draw("rdc-r:z:r","","colz");
170   c->SaveAs(Form("%s_rRes.png",addToName.Data()));
171   //
172   c=GetCanvas(addToName+"_phiRes");
173   t->SetAlias("phiFix","-((phidc-phi)>4)*2*TMath::Pi()+((phidc-phi)<-4)*2*TMath::Pi()");
174   t->Draw("phidc-phi+phiFix:z:r","","colz");
175   c->SaveAs(Form("%s_phiRes.png",addToName.Data()));
176   //
177   c=GetCanvas(addToName+"_rphiRes");
178   t->Draw("(phidc*rdc)-(phi*r):z+(r-84)/(254-84)*18:r","abs(phidc-phi)<1","colz");
179   c->SaveAs(Form("%s_rphiRes.png",addToName.Data()));
180
181   TCanvas *c2=0x0;
182   c2=GetCanvas(addToName+"_Res_1D");
183   c2->Divide(2,2);
184   
185   c2->cd(1);
186   t->Draw("zdc-z","","");
187   //
188   c2->cd(2);
189   t->Draw("rdc-r","","");
190   //
191   c2->cd(3);
192   t->SetAlias("phiFix","-((phidc-phi)>4)*2*TMath::Pi()+((phidc-phi)<-4)*2*TMath::Pi()");
193   t->SetAlias("phiRes","phidc-phi+phiFix");
194   t->Draw("phiRes","","");
195   //
196   c2->cd(4);
197   t->Draw("(phidc*rdc)-(phi*r)","abs(phidc-phi)<1","");
198
199   c2->SaveAs(Form("%s_Res_1D.png",addToName.Data()));
200   
201   f.Close();
202 }
203
204 void makeHistosDist(TString addToName) {
205   TString filename; //("test_");
206   filename.Append(addToName.Data());
207   filename.Append(".root");
208   TFile f(filename.Data());
209   gROOT->cd();
210   TTree *t=(TTree*)f.Get("t");
211   gStyle->SetTitleX(0.18);
212   gStyle->SetTitleW(1-.18-.1);
213   
214   t->SetMarkerStyle(20);
215   t->SetMarkerSize(.8);
216   
217   TCanvas *c=0x0;
218   c=GetCanvas(addToName+"_zResDist");
219   t->Draw("zd-z:z:r","","colz");
220   c->SaveAs(Form("%s_zResDist.png",addToName.Data()));
221   //
222   c=GetCanvas(addToName+"_rResDist");
223   t->Draw("rd-r:z:r","","colz");
224   c->SaveAs(Form("%s_rResDist.png",addToName.Data()));
225   //
226   c=GetCanvas(addToName+"_phiResDist");
227   t->Draw("phid-phi:z:r","abs(phid-phi)<1","colz");
228   c->SaveAs(Form("%s_phiResDist.png",addToName.Data()));
229   //
230   c=GetCanvas(addToName+"_rphiResDist");
231   t->Draw("(phid*rd)-(phi*r):z+(r-84)/(254-84)*18:r","abs(phid-phi)<1","colz");
232   c->SaveAs(Form("%s_rphiResDist.png",addToName.Data()));
233   
234   f.Close();
235 }
236
237
238 void makeAllHistos() {
239   makeHistos("LUT_05");
240   makeHistos("LUT_10");
241   makeHistos("LUT_20");
242   makeHistos("LUT_25");
243   makeHistos("LUT_30");
244   makeHistos("LUT_35");
245   makeHistos("LUT_40");
246   
247 }
248
249 void makeAllHistosDist() {
250   makeHistosDist("LUT_05");
251   makeHistosDist("LUT_10");
252   makeHistosDist("LUT_20");
253   makeHistosDist("LUT_25");
254   makeHistosDist("LUT_30");
255   makeHistosDist("LUT_35");
256   makeHistosDist("LUT_40");
257   
258 }
259
260 TCanvas *GetCanvas(TString addToName)
261 {
262   TString cName(addToName);
263   cName.Prepend("c_");
264   TCanvas *c=(TCanvas*)gROOT->GetListOfCanvases()->FindObject(cName.Data());
265   if (!c) c=new TCanvas(cName.Data(),addToName.Data());
266   c->Clear();
267   c->cd();
268   return c;
269 }
270
271
272