1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 // This implementation AliTPCExB is using an aproximate calculation of the ExB
18 // effect. Therefore the drift ODE is Taylor expanded and only the first
19 // order part is taken.
22 // The ExB correction map is stored in the calib DB
23 // To test current version:
25 char *storage = "local://OCDBres"
27 AliCDBManager::Instance()->SetDefaultStorage(storage);
28 AliCDBManager::Instance()->SetRun(RunNumber)
29 AliTPCExBFirst * exb = AliTPCcalibDB::Instance()->GetExB();
31 // exb->TestExB("exb.root");
32 // TFile f("exb.root");
33 //positions->Draw("drphi");
39 //#include "AliFieldMap.h"
41 #include "TTreeStream.h"
42 #include "AliTPCExBFirst.h"
44 ClassImp(AliTPCExBFirst)
46 const Double_t AliTPCExBFirst::fgkEM=1.602176487e-19/9.10938215e-31;
47 const Double_t AliTPCExBFirst::fgkDriftField=-40.e3;
49 AliTPCExBFirst::AliTPCExBFirst()
51 fkNX(0),fkNY(0),fkNZ(0),
52 fkXMin(-250.),fkXMax(250.),fkYMin(-250.),fkYMax(250.),
53 fkZMin(-250.),fkZMax(250.),
55 fkMeanBx(0),fkMeanBy(0),fkMeanBz(0.) {
62 AliTPCExBFirst::AliTPCExBFirst(const AliMagF *bField,
63 Double_t driftVelocity,
64 Int_t nx,Int_t ny,Int_t nz)
65 : fDriftVelocity(driftVelocity),
66 fkNX(nx),fkNY(ny),fkNZ(nz),
67 fkXMin(-250.),fkXMax(250.),fkYMin(-250.),fkYMax(250.),
68 fkZMin(-250.),fkZMax(250.),
70 fkMeanBx(0),fkMeanBy(0),fkMeanBz(0.) {
72 // The constructor. One has to supply a magnetic field and an (initial)
73 // drift velocity. Since some kind of lookuptable is created the
74 // number of its meshpoints can be supplied.
76 // ConstructCommon(0,bField);
77 ConstructCommon(bField);
82 AliTPCExBFirst::AliTPCExBFirst(const AliFieldMap *bFieldMap,
83 Double_t driftVelocity)
84 : fDriftVelocity(driftVelocity),
85 fkNX(0),fkNY(0),fkNZ(0),
86 fkXMin(-250.),fkXMax(250.),fkYMin(-250.),fkYMax(250.),
87 fkZMin(-250.),fkZMax(250.),
89 fkMeanBx(0),fkMeanBy(0),fkMeanBz(0.) {
91 // The constructor. One has to supply a field map and an (initial)
95 fkXMin=bFieldMap->Xmin()
96 -TMath::Ceil( (bFieldMap->Xmin()+250.0)/bFieldMap->DelX())
98 fkXMax=bFieldMap->Xmax()
99 -TMath::Floor((bFieldMap->Xmax()-250.0)/bFieldMap->DelX())
101 fkYMin=bFieldMap->Ymin()
102 -TMath::Ceil( (bFieldMap->Ymin()+250.0)/bFieldMap->DelY())
104 fkYMax=bFieldMap->Ymax()
105 -TMath::Floor((bFieldMap->Ymax()-250.0)/bFieldMap->DelY())
107 fkZMin=bFieldMap->Zmin()
108 -TMath::Ceil( (bFieldMap->Zmin()+250.0)/bFieldMap->DelZ())
110 fkZMax=bFieldMap->Zmax()
111 -TMath::Floor((bFieldMap->Zmax()-250.0)/bFieldMap->DelZ())
114 fkNX=static_cast<Int_t>((fkXMax-fkXMin)/bFieldMap->DelX()+1.1);
115 fkNY=static_cast<Int_t>((fkYMax-fkYMin)/bFieldMap->DelY()+1.1);
116 fkNZ=static_cast<Int_t>((fkZMax-fkZMin)/bFieldMap->DelZ()+1.1);
118 ConstructCommon(bFieldMap,0);
122 AliTPCExBFirst::~AliTPCExBFirst() {
124 // destruct the poor object.
130 void AliTPCExBFirst::Correct(const Double_t *position,Double_t *corrected) {
132 // correct for the distortion
135 GetMeanFields(position[0],position[1],position[2],&bx,&by);
136 if (position[2]>0.) {
138 GetMeanFields(position[0],position[1],250.,&bxe,&bye);
139 if (position[2]!=250.) {
140 bx=(500.*bxe-(position[2]+250.)*bx)/(250.-position[2]);
141 by=(500.*bye-(position[2]+250.)*by)/(250.-position[2]);
149 Double_t mu=fDriftVelocity/fgkDriftField;
150 Double_t wt=mu*fkMeanBz;
152 corrected[0]=mu*(wt*bx-by)/(1.+wt*wt);
153 corrected[1]=mu*(wt*by+bx)/(1.+wt*wt);
155 if (position[2]>0.) {
156 corrected[0]*=(250.-position[2]);
157 corrected[1]*=(250.-position[2]);
160 corrected[0]*=(-250.-position[2]);
161 corrected[1]*=(-250.-position[2]);
164 corrected[0]=position[0]-corrected[0];
165 corrected[1]=position[1]-corrected[1];
166 corrected[2]=position[2];
169 void AliTPCExBFirst::TestThisBeautifulObject(const char* fileName) {
171 // well, as the name sais...
173 TTreeSRedirector ts(fileName);
175 for (x[0]=-250.;x[0]<=250.;x[0]+=10.)
176 for (x[1]=-250.;x[1]<=250.;x[1]+=10.)
177 for (x[2]=-250.;x[2]<=250.;x[2]+=10.) {
180 Double_t r=TMath::Sqrt(x[0]*x[0]+x[1]*x[1]);
181 Double_t rd=TMath::Sqrt(d[0]*d[0]+d[1]*d[1]);
183 Double_t phi=TMath::ATan2(x[0],x[1]);
184 Double_t phid=TMath::ATan2(d[0],d[1]);
185 Double_t dphi=phi-phid;
186 if (dphi<0.) dphi+=TMath::TwoPi();
187 if (dphi>TMath::Pi()) dphi=TMath::TwoPi()-dphi;
188 Double_t drphi=r*dphi;
189 Double_t dx=x[0]-d[0];
190 Double_t dy=x[1]-d[1];
191 Double_t dz=x[2]-d[2];
208 void AliTPCExBFirst::ConstructCommon(//const AliFieldMap *bFieldMap,
209 const AliMagF *bField) {
211 // THIS IS PRIVATE! (a helper for the constructor)
213 fkNMean=fkNX*fkNY*fkNZ;
214 fkMeanBx=new Double_t[fkNMean];
215 fkMeanBy=new Double_t[fkNMean];
220 for (int i=0;i<fkNX;++i) {
221 x[0]=fkXMin+i*(fkXMax-fkXMin)/(fkNX-1);
222 for (int j=0;j<fkNY;++j) {
223 x[1]=fkYMin+j*(fkYMax-fkYMin)/(fkNY-1);
224 Double_t r=TMath::Sqrt(x[0]*x[0]+x[1]*x[1]);
225 Double_t bx=0.,by=0.;
226 for (int k=0;k<fkNZ;++k) {
227 x[2]=fkZMin+k*(fkZMax-fkZMin)/(fkNZ-1);
229 // the x is not const in the Field function...
231 for (int l=0;l<3;++l) xt[l]=x[l];
232 // that happens due to the lack of a sophisticated class design:
234 // bFieldMap->Field(xt,b);
236 ((AliMagF*)bField)->Field(xt,b);
239 fkMeanBx[(k*fkNY+j)*fkNX+i]=bx;
240 fkMeanBy[(k*fkNY+j)*fkNX+i]=by;
241 if (90.<=r&&r<=250.) {
252 void AliTPCExBFirst::GetMeanFields(Double_t rx,Double_t ry,Double_t rz,
253 Double_t *Bx,Double_t *By) const {
255 // THIS IS PRIVATE! (calculates the mean field utilising a lookup table)
257 Double_t x=(fkNX-1)*(rx-fkXMin)/(fkXMax-fkXMin);
258 Int_t xi1=static_cast<Int_t>(x);
259 xi1=TMath::Max(TMath::Min(xi1,fkNX-2),0);
262 Double_t dx1=(xi2-x);
264 Double_t y=(fkNY-1)*(ry-fkYMin)/(fkYMax-fkYMin);
265 Int_t yi1=static_cast<Int_t>(y);
266 yi1=TMath::Max(TMath::Min(yi1,fkNY-2),0);
269 Double_t dy1=(yi2-y);
271 Double_t z=(fkNZ-1)*(rz-fkZMin)/(fkZMax-fkZMin);
272 Int_t zi1=static_cast<Int_t>(z);
273 zi1=TMath::Max(TMath::Min(zi1,fkNZ-2),0);
277 double s0x=fkMeanBx[yi1*fkNX+xi1]*dx1*dy1
278 +fkMeanBx[yi2*fkNX+xi1]*dx1*dy
279 +fkMeanBx[yi1*fkNX+xi2]*dx *dy1
280 +fkMeanBx[yi2*fkNX+xi2]*dx *dy;
281 double s0y=fkMeanBy[yi1*fkNX+xi1]*dx1*dy1
282 +fkMeanBy[yi2*fkNX+xi1]*dx1*dy
283 +fkMeanBy[yi1*fkNX+xi2]*dx *dy1
284 +fkMeanBy[yi2*fkNX+xi2]*dx *dy;
288 snmx=fkMeanBx[(zi0*fkNY+yi1)*fkNX+xi1]*dx1*dy1
289 +fkMeanBx[(zi0*fkNY+yi2)*fkNX+xi1]*dx1*dy
290 +fkMeanBx[(zi0*fkNY+yi1)*fkNX+xi2]*dx *dy1
291 +fkMeanBx[(zi0*fkNY+yi2)*fkNX+xi2]*dx *dy;
292 snmy=fkMeanBy[(zi0*fkNY+yi1)*fkNX+xi1]*dx1*dy1
293 +fkMeanBy[(zi0*fkNY+yi2)*fkNX+xi1]*dx1*dy
294 +fkMeanBy[(zi0*fkNY+yi1)*fkNX+xi2]*dx *dy1
295 +fkMeanBy[(zi0*fkNY+yi2)*fkNX+xi2]*dx *dy;
299 double snx=fkMeanBx[(zi1*fkNY+yi1)*fkNX+xi1]*dx1*dy1
300 +fkMeanBx[(zi1*fkNY+yi2)*fkNX+xi1]*dx1*dy
301 +fkMeanBx[(zi1*fkNY+yi1)*fkNX+xi2]*dx *dy1
302 +fkMeanBx[(zi1*fkNY+yi2)*fkNX+xi2]*dx *dy;
303 double sny=fkMeanBy[(zi1*fkNY+yi1)*fkNX+xi1]*dx1*dy1
304 +fkMeanBy[(zi1*fkNY+yi2)*fkNX+xi1]*dx1*dy
305 +fkMeanBy[(zi1*fkNY+yi1)*fkNX+xi2]*dx *dy1
306 +fkMeanBy[(zi1*fkNY+yi2)*fkNX+xi2]*dx *dy;
307 double snpx=fkMeanBx[(zi2*fkNY+yi1)*fkNX+xi1]*dx1*dy1
308 +fkMeanBx[(zi2*fkNY+yi2)*fkNX+xi1]*dx1*dy
309 +fkMeanBx[(zi2*fkNY+yi1)*fkNX+xi2]*dx *dy1
310 +fkMeanBx[(zi2*fkNY+yi2)*fkNX+xi2]*dx *dy;
311 double snpy=fkMeanBy[(zi2*fkNY+yi1)*fkNX+xi1]*dx1*dy1
312 +fkMeanBy[(zi2*fkNY+yi2)*fkNX+xi1]*dx1*dy
313 +fkMeanBy[(zi2*fkNY+yi1)*fkNX+xi2]*dx *dy1
314 +fkMeanBy[(zi2*fkNY+yi2)*fkNX+xi2]*dx *dy;
318 *Bx=0.5*(((snpx-2.*snx+snmx)*dz+2.*(snx-snmx))*dz+snx-s0x+snmx);
319 *By=0.5*(((snpy-2.*sny+snmy)*dz+2.*(sny-snmy))*dz+sny-s0y+snmy);
320 //TODO: make this nice
321 if (TMath::Abs(z)>0.001) {