6 #include "TTreeStream.h"
10 #include "AliTPCCorrection.h"
11 #include "AliTPCCorrectionLookupTable.h"
13 void makeComparisonTree(TString filename, TString addToName)
15 TFile fn(filename.Data());
17 AliTPCCorrectionLookupTable *fTPCCorrection=(AliTPCCorrectionLookupTable*)fn.Get("map");
19 // fTPCCorrection->BuildExactInverse();
21 // TFile f("/tmp/corrTest.Root","recreate");
22 // fTPCCorrection->Write("map");
25 TString outFile=addToName;
26 outFile.Append(".root");
27 TTreeSRedirector *sred=new TTreeSRedirector(outFile.Data());
29 Float_t dx[3]={0,0,0};
31 for (Float_t iz=-245; iz<=245; iz+=10) {
32 Short_t roc=(iz>=0)?0:18;
33 for (Float_t ir=86; ir<250; ir+=10) {
34 for (Float_t iphi=0; iphi<TMath::TwoPi(); iphi+=10*TMath::DegToRad()){
35 Float_t x=ir*(Float_t)TMath::Cos(iphi);
36 Float_t y=ir*(Float_t)TMath::Sin(iphi);
37 Float_t x3[3]={x,y,iz};
46 fTPCCorrection->GetDistortion(x3,roc,dx);
49 Float_t zd = iz+dx[2];
50 Float_t rd = TMath::Sqrt(xd*xd+yd*yd);
51 Float_t phid = TMath::ATan2(yd,xd);
52 if (phid<0) phid+=TMath::TwoPi();
60 // correct back distorted point
61 Float_t xd3[3]={xd,yd,zd};
63 fTPCCorrection->GetCorrection(xd3,roc,dx);
64 Float_t xdc = xd+dx[0];
65 Float_t ydc = yd+dx[1];
66 Float_t zdc = zd+dx[2];
67 Float_t rdc = TMath::Sqrt(xdc*xdc+ydc*ydc);
68 Float_t phidc = TMath::ATan2(ydc,xdc);
69 if (phidc<0) phidc+=TMath::TwoPi();
78 // write current point
89 void makeAllComparisonTrees()
91 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2_eps5_50kHz_precal.lookup.root","LUT_05");
92 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2_eps10_50kHz_precal.lookup.root","LUT_10");
93 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2_eps20_50kHz_precal.lookup.root","LUT_20");
96 void makeAllComparisonTreesNew()
98 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps5_50kHz_precal.lookup.root","LUT_05");
99 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps10_50kHz_precal.lookup.root","LUT_10");
100 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_20");
101 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_25");
102 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_30");
103 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_35");
104 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_40");
107 void makeAllComparisonTreesOld()
109 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps5_50kHz_precal.lookup.root","LUT_05");
110 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps10_50kHz_precal.lookup.root","LUT_10");
111 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_20");
112 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_25");
113 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_30");
114 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_35");
115 makeComparisonTree("$ALICE_ROOT/TPC/Calib/maps/old/SC_NeCO2N2_eps20_50kHz_precal.lookup.root","LUT_40");
118 TCanvas *GetCanvas(TString addToName);
120 void makeHistos(TString addToName) {
121 TString fileName; //("test_");
122 fileName.Append(addToName.Data());
123 fileName.Append(".root");
124 TFile f(fileName.Data());
126 TTree *t=(TTree*)f.Get("t");
127 gStyle->SetTitleX(0.18);
128 gStyle->SetTitleW(1-.18-.1);
130 t->SetMarkerStyle(20);
131 t->SetMarkerSize(.8);
134 c=GetCanvas(addToName+"_zRes");
135 t->Draw("zdc-z:z:r","","colz");
136 c->SaveAs(Form("%s_zRes.png",addToName.Data()));
138 c=GetCanvas(addToName+"_rRes");
139 t->Draw("rdc-r:z:r","","colz");
140 c->SaveAs(Form("%s_rRes.png",addToName.Data()));
142 c=GetCanvas(addToName+"_phiRes");
143 t->SetAlias("phiFix","-((phidc-phi)>4)*2*TMath::Pi()+((phidc-phi)<-4)*2*TMath::Pi()");
144 t->Draw("phidc-phi+phiFix:z:r","","colz");
145 c->SaveAs(Form("%s_phiRes.png",addToName.Data()));
147 c=GetCanvas(addToName+"_rphiRes");
148 t->Draw("(phidc*rdc)-(phi*r):z+(r-84)/(254-84)*18:r","abs(phidc-phi)<1","colz");
149 c->SaveAs(Form("%s_rphiRes.png",addToName.Data()));
152 c2=GetCanvas(addToName+"_Res_1D");
156 t->Draw("zdc-z","","");
159 t->Draw("rdc-r","","");
162 t->SetAlias("phiFix","-((phidc-phi)>4)*2*TMath::Pi()+((phidc-phi)<-4)*2*TMath::Pi()");
163 t->Draw("phidc-phi+phiFix","","");
166 t->Draw("(phidc*rdc)-(phi*r)","abs(phidc-phi)<1","");
168 c2->SaveAs(Form("%s_Res_1D.png",addToName.Data()));
173 void makeHistosDist(TString addToName) {
174 TString fileName; //("test_");
175 fileName.Append(addToName.Data());
176 fileName.Append(".root");
177 TFile f(fileName.Data());
179 TTree *t=(TTree*)f.Get("t");
180 gStyle->SetTitleX(0.18);
181 gStyle->SetTitleW(1-.18-.1);
183 t->SetMarkerStyle(20);
184 t->SetMarkerSize(.8);
187 c=GetCanvas(addToName+"_zResDist");
188 t->Draw("zd-z:z:r","","colz");
189 c->SaveAs(Form("%s_zResDist.png",addToName.Data()));
191 c=GetCanvas(addToName+"_rResDist");
192 t->Draw("rd-r:z:r","","colz");
193 c->SaveAs(Form("%s_rResDist.png",addToName.Data()));
195 c=GetCanvas(addToName+"_phiResDist");
196 t->Draw("phid-phi:z:r","abs(phid-phi)<1","colz");
197 c->SaveAs(Form("%s_phiResDist.png",addToName.Data()));
199 c=GetCanvas(addToName+"_rphiResDist");
200 t->Draw("(phid*rd)-(phi*r):z+(r-84)/(254-84)*18:r","abs(phid-phi)<1","colz");
201 c->SaveAs(Form("%s_rphiResDist.png",addToName.Data()));
207 void makeAllHistos() {
208 makeHistos("LUT_05");
209 makeHistos("LUT_10");
210 makeHistos("LUT_20");
211 makeHistos("LUT_25");
212 makeHistos("LUT_30");
213 makeHistos("LUT_35");
214 makeHistos("LUT_40");
218 void makeAllHistosDist() {
219 makeHistosDist("LUT_05");
220 makeHistosDist("LUT_10");
221 makeHistosDist("LUT_20");
222 makeHistosDist("LUT_25");
223 makeHistosDist("LUT_30");
224 makeHistosDist("LUT_35");
225 makeHistosDist("LUT_40");
229 TCanvas *GetCanvas(TString addToName)
231 TString cName(addToName);
233 TCanvas *c=(TCanvas*)gROOT->GetListOfCanvases()->FindObject(cName.Data());
234 if (!c) c=new TCanvas(cName.Data(),addToName.Data());