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