]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/CalibMacros/AliTPCCorrectionDemo.C
Adding macro for ExB global fit
[u/mrichter/AliRoot.git] / TPC / CalibMacros / AliTPCCorrectionDemo.C
1 /*
2  gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/ -I$ALICE_ROOT/include -I$ALICE_ROOT/STEER -I$ALICE_ROOT/TPC -I$ALICE_ROOT/ITS -I$ALICE_ROOT/TRD -I$ALICE_ROOT/TOF -I$ALICE_ROOT/RAW -I$ALICE_ROOT/STAT");
3
4  */
5 #if !defined(__CINT__) || defined(__MAKECINT__)
6 #include "THnSparse.h"
7 #include "TLatex.h"
8 #include "TCanvas.h"
9 #include "TLegend.h"
10 #include "TSystem.h"
11 #include "TFile.h"
12 #include "TChain.h"
13 #include "TCut.h"
14 #include "TH3.h"
15 #include "TH2F.h"
16 #include "TProfile3D.h"
17 #include "TMath.h" 
18 #include "TVectorD.h"
19 #include "TMatrixD.h"
20 #include "TStatToolkit.h"
21 #include "TTreeStream.h"
22 #include "AliExternalTrackParam.h"
23 #include "AliESDfriend.h"
24 #include "AliTPCcalibTime.h"
25 #include "TROOT.h"
26 #include "AliXRDPROOFtoolkit.h"
27 #include "AliTPCCorrection.h"
28 #include "AliTPCExBTwist.h"
29 #include "AliTPCGGVoltError.h"
30 #include "AliTPCComposedCorrection.h"
31 #include "AliTPCExBConical.h"
32 #include "TPostScript.h"
33 #include "TStyle.h"
34 #include "AliTrackerBase.h"
35 #include "TGraph.h"
36 #include "AliCDBManager.h"
37 #include "AliTPCExBBShape.h"
38 #include "TRandom.h"
39 #include "AliGeomManager.h"
40
41 #endif
42
43
44 void AliTPCCorrectionDemo() {
45
46   //
47   // This is a Demo function of the general class AliTPCCorrection, which is used for 
48   // general space point correction due to different effects.
49   // The effects used in this Demo are:
50   //   1. ExB twist - general offset of the TPC axis in comparison to the B field axis
51   //   2. GG error (Gating Grid volt. error) - not perfectly aligned GG voltage (in terms of voltage)
52   //   3. ExBBShape - B field shape correction of the secound order
53   //
54   // See class descriptions for further details 
55   //
56   // Authors: Magnus Mager, Stefan Rossegger, Jim Thomas
57   //
58   //
59   // omegaTau (wt) of the langevin equation
60   // This is a function of the drift vel., the magnetic and electric field
61   // e.g. vd=2.6 cm/usc; Ez=400 V/cm; Bz=0.5 T
62   // wt =  -10.0*(Bz*10)*vd/Ez = -0.325 
63
64   Double_t vdrift = 2.6; // [cm/us]   // to be updated: per second (ideally)
65   Double_t bzField = -0.5; // [Tesla] // to be updated: per run
66   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
67   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
68
69   // Correction Terms for effective omegaTau; obtained by a laser calibration run
70   Double_t T1 = 0.9;
71   Double_t T2 = 1.5;
72
73   AliMagF mag("mag","mag");
74
75   AliTPCExBTwist twist;
76   twist.SetXTwist(0.001);
77   
78   AliTPCGGVoltError GGerror;
79   GGerror.SetDeltaVGGA(50.);
80   GGerror.SetDeltaVGGC(50.);
81   GGerror.InitGGVoltErrorDistortion();
82
83   AliTPCExBBShape exb;
84   exb.SetBField(&mag);
85
86   TObjArray cs;
87   cs.Add(&twist);
88   cs.Add(&GGerror);
89   cs.Add(&exb);
90
91   AliTPCComposedCorrection cc;
92   cc.SetCorrections(&cs);
93   cc.SetOmegaTauT1T2(wt,T1,T2);
94   //cc.SetMode(1);
95
96   cc.Print("DA"); // Print used correction classes
97
98   TCanvas *c=new TCanvas;  // Plots
99   c->Divide(2,2);
100   c->cd(1);twist.CreateHistoDRPhiinZR(1.,100,100)->Draw("surf2");
101   c->cd(2);GGerror.CreateHistoDRPhiinZR(1.,100,100)->Draw("surf2");
102   c->cd(3);exb.CreateHistoDRPhiinZR(1.)->Draw("surf2");
103   c->cd(4);cc.CreateHistoDRPhiinZR(1.)->Draw("surf2");
104 }
105
106
107
108 void MakeDistortionMap(){
109   //
110   // make distortiona map example for specific transformation 
111   //
112   Int_t run=0;
113   TTreeSRedirector *pcstream =  new TTreeSRedirector("distort.root");
114   Double_t vdrift = 2.6; // [cm/us]   // to be updated: per second (ideally)
115   Double_t bzField = -0.5; // [Tesla] // to be updated: per run
116   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
117   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
118
119   // Correction Terms for effective omegaTau; obtained by a laser calibration run
120   Double_t T1 = 0.9;
121   Double_t T2 = 1.5;
122   AliMagF* magF= new AliMagF("Maps","Maps", bzField/0.5, 1.);
123   TGeoGlobalMagField::Instance()->SetField(magF);
124
125   AliCDBManager::Instance()->SetDefaultStorage("local:////$ALICE_ROOT/OCDB");
126   AliCDBManager::Instance()->SetRun(run);
127   AliGeomManager::LoadGeometry();
128   AliGeomManager::ApplyAlignObjsFromCDB("GRP ITS TPC");
129
130
131   AliTPCExBTwist twist;
132   twist.SetName("twistX001");
133   twist.SetXTwist(0.001);
134   
135   AliTPCGGVoltError GGerror;
136   GGerror.SetDeltaVGGA(50.);
137   GGerror.SetDeltaVGGC(50.);
138   GGerror.InitGGVoltErrorDistortion();
139   GGerror.SetName("ggoffsetA50C50");
140   AliTPCExBBShape exb;
141   exb.SetBField(magF);
142   exb.SetName("ExB");
143   TObjArray cs;
144   cs.Add(&twist);
145   cs.Add(&GGerror);
146   cs.Add(&exb);
147   
148   AliTPCComposedCorrection *cc = new AliTPCComposedCorrection;
149   cc->SetCorrections(&cs);
150   cc->SetOmegaTauT1T2(wt,T1,T2);
151   //cc->SetMode(1);
152   cc->SetName("composed");
153   Double_t par[5]={0,0,0,0,0};
154   Double_t cov[15]={0,0,0,0,0,0,0,0,0,0,0,0,0,0};
155   for (Double_t theta=-1.; theta<1; theta+=0.1){
156     for (Double_t alpha=-3.14; alpha<3.14; alpha+=0.2){
157       par[0]=0;
158       par[1]=theta*85;
159       par[2]=0;
160       par[3]=theta;
161       par[4]=(gRandom->Rndm()-0.5);
162       AliExternalTrackParam trackIn(85,alpha,par,cov);
163       //FitDistortedTrack(&trackIn, cc, 85, -1,pcstream);
164       twist.FitDistortedTrack(trackIn, 85, -1,pcstream);
165       GGerror.FitDistortedTrack(trackIn, 85, -1,pcstream);
166       exb.FitDistortedTrack(trackIn, 85, -1,pcstream);
167       cc->FitDistortedTrack(trackIn, 85, -1,pcstream);      
168     }
169   }
170   delete pcstream;
171 }
172
173 void DrawDistortionMap(){
174   //
175   // Example -drawing of distortion maps
176   //
177   TFile f("distort.root");
178   TTree *   fitDistorttwistX001= (TTree*)f.Get("fitDistorttwistX001");
179   TTree *   fitDistortggoffsetA50C50=(TTree*)f.Get("fitDistortggoffsetA50C50");
180   TTree *   fitDistortExB=(TTree*)f.Get("fitDistortExB");
181   TTree *   fitDistortcomposed=(TTree*)f.Get("fitDistortcomposed");
182   fitDistortExB->SetMarkerColor(1);
183   fitDistorttwistX001->SetMarkerColor(2);
184   fitDistortggoffsetA50C50->SetMarkerColor(4);
185   fitDistortExB->SetMarkerStyle(20);
186   fitDistorttwistX001->SetMarkerStyle(21);
187   fitDistortggoffsetA50C50->SetMarkerStyle(22);
188   //
189   // example draw delta local y/r-phi as function of fi
190   Int_t entries=0;
191
192   entries = fitDistortExB->Draw("track1.fP[0]-track0.fP[0]:track0.fAlpha","","");
193   TGraph * grY0= new TGraph(entries, fitDistortExB->GetV2(),fitDistortExB->GetV1());
194   entries = fitDistortggoffsetA50C50->Draw("track1.fP[0]-track0.fP[0]:track0.fAlpha","","");
195   TGraph * grY1= new TGraph(entries, fitDistortggoffsetA50C50->GetV2(),fitDistortggoffsetA50C50->GetV1());
196   entries = fitDistorttwistX001->Draw("track1.fP[0]-track0.fP[0]:track0.fAlpha","","");
197   TGraph * grY2= new TGraph(entries, fitDistorttwistX001->GetV2(),fitDistorttwistX001->GetV1());
198   entries = fitDistortcomposed->Draw("track1.fP[0]-track0.fP[0]:track0.fAlpha","","");
199   TGraph * grY3= new TGraph(entries, fitDistortcomposed->GetV2(),fitDistortcomposed->GetV1());
200   grY0->SetMinimum(-0.4);
201   grY0->SetMaximum(0.4);
202   grY0->SetMarkerStyle(20), grY0->SetMarkerColor(1);
203   grY1->SetMarkerStyle(21), grY1->SetMarkerColor(2);
204   grY2->SetMarkerStyle(22), grY2->SetMarkerColor(4);
205   grY3->SetMarkerStyle(24), grY3->SetMarkerColor(5); 
206   grY0->Draw("ap");
207   grY1->Draw("p");
208   grY2->Draw("p");
209   grY3->Draw("p");
210 }
211