]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCExBEffectiveSector.cxx
Add macro to create OCDB entry for ion tail
[u/mrichter/AliRoot.git] / TPC / AliTPCExBEffectiveSector.cxx
CommitLineData
40787a37 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16////////////////////////////////////////////////////////////////////////////
17// //
18// AliTPCExBEffectiveSector class //
19// Correct for the rest of ExB effect which are not covered yet by physical models
20//
21// Motivation:
22// ExB correction:
23// dr = c0* integral(Er/Ez) + c1* integral(Erphi/Ez)
24// drphi = -c1* integral(Er/Ez) + c0* integral(Erphi/Ez)
25// Where:
26// wt = Bz*(k*vdrift/E) ~ 0.3 at B=0.5 T
27// c0 = 1/(1+T2*T2*wt*wt)
28// c1 = T1*wt/(1+T1*T1*wt*wt)
29//
30//
31// 3 correction maps 0 implemented as histogram used
32// R-Phi correction map obtained minimizing residuals betwee the track
33// and space points (AliTPCcalibAlign class). Track is defined using
34// the points from the refernce plain at the middle of the TPC
35// and vertex
36// Corrected primar tracks straight pointing to the primary vertex
37//
38// R distortion - obtained using the cluster residuals in the setup with
39// plus and minus field
40// Only high momenta tracks used for this calibration (1 GeV threshold)
41// drphi_plus-drphi_minus=-2*c1 integral(Er/Ez)
42// - Erphi/Ez cancels
43////////////////////////////////////////////////////////////////////////////
44#include "AliMagF.h"
45#include "TGeoGlobalMagField.h"
46#include "AliTPCcalibDB.h"
47#include "AliTPCParam.h"
48#include "AliLog.h"
49
50#include "TMath.h"
51#include "AliTPCROC.h"
52#include "TFile.h"
53#include "TAxis.h"
54#include "TTree.h"
55#include "TTreeStream.h"
56#include "THnSparse.h"
7f00a31e 57#include "THnBase.h"
40787a37 58#include "TProfile.h"
59#include "TH2F.h"
60#include "TH3F.h"
61#include "TROOT.h"
62#include "AliTPCExBEffectiveSector.h"
63ClassImp(AliTPCExBEffectiveSector)
64
65AliTPCExBEffectiveSector::AliTPCExBEffectiveSector()
66 : AliTPCCorrection("ExB_effectiveSector","ExB effective sector"),
67 fC0(1.),fC1(0.),
68 fCorrectionR(0), // radial correction
69 fCorrectionRPhi(0), // r-phi correction
70 fCorrectionZ(0) // z correction
71{
72 //
73 // default constructor
74 //
75}
76
77AliTPCExBEffectiveSector::~AliTPCExBEffectiveSector() {
78 //
79 // default destructor
80 //
81 delete fCorrectionR; // radial correction
82 delete fCorrectionRPhi; // r-phi correction
83 delete fCorrectionZ; // z correction
84}
85
86
87
88void AliTPCExBEffectiveSector::Init() {
89 //
90 // Initialization funtion
91 //
92
93 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
94 if (!magF) AliError("Magneticd field - not initialized");
95 Double_t bzField = magF->SolenoidField()/10.; //field in T
96 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
97 if (!param) AliError("Parameters - not initialized");
98 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
99 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
100 Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
101 // Correction Terms for effective omegaTau; obtained by a laser calibration run
102 SetOmegaTauT1T2(wt,fT1,fT2);
103}
104
105void AliTPCExBEffectiveSector::Update(const TTimeStamp &/*timeStamp*/) {
106 //
107 // Update function
108 //
109 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
110 if (!magF) AliError("Magneticd field - not initialized");
111 Double_t bzField = magF->SolenoidField()/10.; //field in T
112 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
113 if (!param) AliError("Parameters - not initialized");
114 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
115 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
116 Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
117 // Correction Terms for effective omegaTau; obtained by a laser calibration run
118 SetOmegaTauT1T2(wt,fT1,fT2);
119}
120
121
122
123void AliTPCExBEffectiveSector::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
124 //
125 // Calculates the correction using the lookup table (histogram) of distortion
126 // The histogram is created as poscl - postrack
127 //
28178a25 128 dx[0]=0;
129 dx[1]=0;
130 dx[2]=0;
131 if (!fCorrectionRPhi) return;
40787a37 132 Double_t phi = TMath::ATan2(x[1],x[0]);
133 Double_t r = TMath::Sqrt(x[1]*x[1]+x[0]*x[0]);
134 Double_t sector = 9.*phi/TMath::Pi();
135 if (sector<0) sector+=18.;
136 Double_t kZ=x[2]/r;
137 //
138 if (kZ>1.2) kZ= 1.2;
139 if (kZ<-1.2) kZ= -1.2;
140 if (roc%36<18) kZ= TMath::Abs(kZ);
141 if (roc%36>=18) kZ=-TMath::Abs(kZ);
142 if (TMath::Abs(kZ)<0.15){
143 kZ = (roc%36<18) ? 0.15:-0.15;
28178a25 144 }
40787a37 145 //
146 Double_t dlR=0;
147 Double_t dlRPhi=0;
148 Double_t dlZ=0;
28178a25 149 Double_t rr=TMath::Max(r,fCorrectionRPhi->GetYaxis()->GetXmin()+0.01);
150 rr=TMath::Min(rr,fCorrectionRPhi->GetYaxis()->GetXmax()-0.01);
151 Double_t kZZ=TMath::Max(kZ,fCorrectionRPhi->GetZaxis()->GetXmin()+0.001);
152 kZZ=TMath::Min(kZZ,fCorrectionRPhi->GetZaxis()->GetXmax()-0.001);
153
40787a37 154 if (fCorrectionRPhi) {
28178a25 155 // dlRPhi= -fCorrectionRPhi->Interpolate(sector,rr,kZZ);
156 dlRPhi= -fCorrectionRPhi->GetBinContent(fCorrectionRPhi->FindBin(sector,rr,kZZ));
40787a37 157 }
158 if (fCorrectionR) {
28178a25 159 // dlR= -fCorrectionR->Interpolate(sector,rr,kZZ);
160 dlR= -fCorrectionR->GetBinContent(fCorrectionR->FindBin(sector,rr,kZZ));
40787a37 161 }
162 if (fCorrectionZ) {
28178a25 163 // dlZ= -fCorrectionZ->Interpolate(sector,rr,kZZ);
164 dlZ= -fCorrectionZ->GetBinContent(fCorrectionZ->FindBin(sector,rr,kZZ));
40787a37 165 }
166 Double_t dr = fC0*dlR + fC1*dlRPhi;
167 Double_t drphi = -fC1*dlR + fC0*dlRPhi;
168 // Calculate distorted position
169 if ( r > 0.0 ) {
170 r = r + dr;
171 phi = phi + drphi/r;
172 }
173 // Calculate correction in cartesian coordinates
174 dx[0] = r * TMath::Cos(phi) - x[0];
175 dx[1] = r * TMath::Sin(phi) - x[1];
176 dx[2] = dlZ;
177
178}
179
180void AliTPCExBEffectiveSector::Print(const Option_t* option) const {
181 //
182 // Print function to check the settings (e.g. the twist in the X direction)
183 // option=="a" prints the C0 and C1 coefficents for calibration purposes
184 //
185
186 TString opt = option; opt.ToLower();
187 printf("%s\t%s\n",GetName(),GetTitle());
188 if (opt.Contains("a")) { // Print all details
189 printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
190 printf(" - C0: %1.4f, C1: %1.4f \n",fC0,fC1);
191 }
192}
193