]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/CalibMacros/AliTPCDistortionFun.C
Update timestamps for new AMANDA simulation (17/02/2015)
[u/mrichter/AliRoot.git] / TPC / CalibMacros / AliTPCDistortionFun.C
CommitLineData
35d81915 1/// \file AliTPCDistortionFun.C
2///
3/// Simple compiled macro for declaration of static distortion function
4/// on top of the AliTPCDistortion class.
5/// Why:
6/// 1. Use static function in the fitting procedure
7/// 2. Usage in TFormual, TF1, Tf2 ... for visualization.
8/// 3. Usage in tree->Draw() for visualization
9/// 4. Simple visualization of fit residuals in multidemension - using tree Draw functionality
10///
11/// Usage:
12/// ~~~
13/// gSystem->AddIncludePath("-I$ALICE_ROOT/TPC");
14/// .L $ALICE_ROOT/TPC/CalibMacros/AliTPCDistortions.cxx+
15/// .L $ALICE_ROOT/TPC/CalibMacros/AliTPCDistortionFun.C+
16/// ~~~
17///
18/// Example:
19///
20/// ~~~{.cpp}
21/// // Draw integrated distortion in local x
22///
23/// TF2 fdistIFCXvZX("fdistIFCXvZX","GetIFCDistortion(y,0,x,0)*sign(x)",-250,250,80,250);
24/// fdistIFCXvZX.SetNpx(200);
25/// fdistIFCXvZX.SetNpy(200);
26/// fdistIFCXvZX.GetXaxis()->SetTitle("Z (cm)");
27/// fdistIFCXvZX.GetYaxis()->SetTitle("local X (cm)");
28/// fdistIFCXvZX->Draw("colz");
29///
30/// // Draw local distortion angle dx/dz in mrad
31/// TF2 fangleIFCXvZX("fangleIFCXvZX","1000*(GetIFCDistortion(y,0,x,0)-GetIFCDistortion(y,0,x-1,0))",-250,250,85,245);
32/// fangleIFCXvZX.SetNpx(200);
33/// fangleIFCXvZX.SetNpy(200);
34/// fangleIFCXvZX.GetXaxis()->SetTitle("Z (cm)");
35/// fangleIFCXvZX.GetYaxis()->SetTitle("local X (cm)");
36/// fangleIFCXvZX->Draw("colz");
37///
38/// TF2 fangleGGXvZX("fangleGGXvZX","1000*(GetGGDistortion(y,0,x,0,1*400,1*400)-GetGGDistortion(y,0,x-1,0,1*400,1*400))*sign(x)",-250,250,85,245);
39/// fangleGGXvZX.SetNpx(200);
40/// fangleGGXvZX.SetNpy(200);
41/// fangleGGXvZX.GetXaxis()->SetTitle("Z (cm)");
42/// fangleGGXvZX.GetYaxis()->SetTitle("local X (cm)");
43/// fangleGGXvZX->Draw("colz");
44/// ~~~
57f936af 45
7f550ed6 46#include "TCanvas.h"
47#include "TF1.h"
48#include "TLegend.h"
57f936af 49#include "AliTPCDistortions.h"
50
7f550ed6 51TObjArray *arrayPic=new TObjArray;
57f936af 52
53Double_t GetIFCDistortion(Double_t lx, Double_t ly, Double_t lz, Int_t icoord, Double_t shift=1.){
35d81915 54 ///
55
57f936af 56 static AliTPCDistortions transform;
7f550ed6 57 static Bool_t doInit=kTRUE;
58 if (doInit){
59 transform.SetIFCShift(1);
60 transform.InitIFCShiftDistortion();
61 doInit=kFALSE;
62 }
63
57f936af 64 Double_t xyzIn[3]={lx, ly, lz};
65 Double_t xyzOut[3]={lx, ly, lz};
57f936af 66 Int_t dummyROC=0;
67 if (lz<0) { dummyROC=36;}
68 transform.UndoIFCShiftDistortion(xyzIn,xyzOut,dummyROC);
69 Double_t result=0;
70 if (icoord<3) result=xyzOut[icoord]-xyzIn[icoord];
71 return result*shift;
72}
7f550ed6 73
74
75
76Double_t GetGGDistortion(Double_t lx, Double_t ly, Double_t lz, Int_t icoord, Double_t deltaVGGA=1., Double_t deltaVGGC=1.){
35d81915 77 /// GG distortion induced distortions
78
7f550ed6 79 static AliTPCDistortions transform;
80 static Bool_t doInit=kTRUE;
81 if (doInit){
82 transform.SetDeltaVGGA(1.);
83 transform.SetDeltaVGGC(1.);
84 transform.InitGGVoltErrorDistortion();
85 doInit=kFALSE;
86 }
87 Double_t xyzIn[3]={lx, ly, lz};
88 Double_t xyzOut[3]={lx, ly, lz};
89 Int_t dummyROC=0;
90 if (lz<0) { dummyROC=36;}
91 transform.UndoGGVoltErrorDistortion(xyzIn,xyzOut,dummyROC);
92 Double_t result=0;
93 if (icoord<3) result=xyzOut[icoord]-xyzIn[icoord];
94 if (lz<0) result*=deltaVGGA;
95 if (lz>0) result*=deltaVGGC;
96 return result;
97}
98
99
100void MakePicIFCDX(){
35d81915 101 ///
102
7f550ed6 103 TCanvas *canvasIFC1D= new TCanvas("IFC radial shift", "IFC radial shift");
104 TLegend *legend = new TLegend(0.1,0.60,0.5,0.9, "Radial distortion due IFC shift 1 mm ");
105
106 for (Int_t iradius=0; iradius<=10; iradius++){
107 Double_t radius=85+iradius*(245-85)/10.;
108 TF1 *f1= new TF1(Form("fIFC_ZX%f",radius),Form("0.1*GetIFCDistortion(%f,0,x,0)*sign(x)",radius),-250,250);
109 f1->SetMaximum( 0.6);
110 f1->SetMinimum(-0.6);
111 f1->SetNpx(200);
112 f1->SetLineColor(1+((20+iradius)%20));
113 f1->SetLineWidth(1);
114 if (iradius==0) f1->Draw("p");
115 f1->Draw("samep");
116 legend->AddEntry(f1,Form("R=%f",radius));
117 }
118 legend->Draw();
119}
120
121void MakePicIFCDXangle(){
35d81915 122 ///
123
7f550ed6 124 TCanvas *canvasIFC1D= new TCanvas("IFC radial shift - angle", "IFC radial shift angle");
125 TLegend *legend = new TLegend(0.1,0.60,0.9,0.9, "Radial distortion due IFC shift 1 mm ");
126
127 for (Int_t iradius=0; iradius<=10; iradius++){
128 Double_t radius=85+iradius*(245-85)/10.;
129 TF1 *f1= new TF1(Form("fIFC_ZX%f",radius),Form("0.1*1000*(GetIFCDistortion(%f,0,x,0)-GetIFCDistortion(%f,0,x+1,0))",radius,radius),-250,250);
130 f1->SetMaximum( 10);
131 f1->SetMinimum(-10);
132 f1->SetNpx(200);
133 f1->SetLineColor(1+(iradius%10));
134 f1->SetLineWidth(1);
135 //f1->SetLineStyle(2+(iradius%6));
136 if (iradius==0) f1->Draw("");
137 f1->Draw("same");
138 legend->AddEntry(f1,Form("R=%f",radius));
139 }
140 legend->Draw();
141}
142