AliRoot Core  edcc906 (edcc906)
AliTPCExB.cxx
Go to the documentation of this file.
1 #include "AliTPCExB.h"
2 #include "TMath.h"
3 //#include "TTreeStream.h"
4 #include "AliMagF.h"
5 #include "TLinearFitter.h"
6 #include "AliTPCcalibDB.h"
7 
8 
27 
29 
31 
32 
33 
35 ClassImp(AliTPCExB)
37 
39  TObject(),
40  fMatBrBz(0), //param matrix Br/Bz
41  fMatBrfiBz(0), //param matrix Br/Bz
42  fMatBrBzI0(0), //param matrix Br/Bz integral z>0
43  fMatBrBzI1(0), //param matrix Br/Bz integral z<0
44  fMatBrfiBzI0(0), //param matrix Br/Bz integral z>0
45  fMatBrfiBzI1(0) //param matrix Br/Bz integral z<0
46 {
47  //
48  // default constructor
49  //
50 }
51 
53  TObject(exb),
54  fMatBrBz(new TVectorD(*(exb.fMatBrBz))), //param matrix Br/Bz
55  fMatBrfiBz(new TVectorD(*(exb.fMatBrfiBz))), //param matrix Br/Bz
56  fMatBrBzI0(new TVectorD(*(exb.fMatBrBzI0))), //param matrix Br/Bz integral z>0
57  fMatBrBzI1(new TVectorD(*(exb.fMatBrBzI1))), //param matrix Br/Bz integral z<0
58  fMatBrfiBzI0(new TVectorD(*(exb.fMatBrfiBzI0))), //param matrix Br/Bz integral z>0
59  fMatBrfiBzI1(new TVectorD(*(exb.fMatBrfiBzI1))) //param matrix Br/Bz integral z<0
60 {
62 
63 }
64 
66 {
68 
69  return *this;
70 }
71 
72 
73 
74 
75 void AliTPCExB::TestExB(const char* fileName) {
78 
79  TTreeSRedirector ts(fileName);
80  Double_t x[3];
81  for (x[0]=-250.;x[0]<=250.;x[0]+=10.)
82  for (x[1]=-250.;x[1]<=250.;x[1]+=10.)
83  for (x[2]=-250.;x[2]<=250.;x[2]+=20.) {
84  Double_t r=TMath::Sqrt(x[0]*x[0]+x[1]*x[1]);
85  if (r<20) continue;
86  if (r>260) continue;
87  Double_t z = x[2];
88  Double_t d[3];
89  Correct(x,d);
90  Double_t rd=TMath::Sqrt(d[0]*d[0]+d[1]*d[1]);
91  Double_t dr=r-rd;
92  Double_t phi=TMath::ATan2(x[1],x[0]);
93  Double_t phid=TMath::ATan2(d[1],d[0]);
94  Double_t dphi=phi-phid;
95  if (dphi<0.) dphi+=TMath::TwoPi();
96  if (dphi>TMath::Pi()) dphi=TMath::TwoPi()-dphi;
97  Double_t drphi=r*dphi;
98  Double_t dx=x[0]-d[0];
99  Double_t dy=x[1]-d[1];
100  Double_t dz=x[2]-d[2];
101  //
102  Double_t bx = GetBx(r,phi,z,0);
103  Double_t by = GetBy(r,phi,z,0);
104  Double_t bz = GetBz(r,phi,z,0);
105  Double_t br = GetBr(r,phi,z,0);
106  Double_t brfi = GetBrfi(r,phi,z,0);
107  //
108  Double_t bxi = GetBxI(r,phi,z,0);
109  Double_t byi = GetByI(r,phi,z,0);
110  Double_t bzi = GetBzI(r,phi,z,0);
111  Double_t bri = GetBrI(r,phi,z,0);
112  Double_t brfii = GetBrfiI(r,phi,z,0);
113 
114  ts<<"positions"<<
115  "x0="<<x[0]<<
116  "x1="<<x[1]<<
117  "x2="<<x[2]<<
118  "dx="<<dx<<
119  "dy="<<dy<<
120  "dz="<<dz<<
121  "r="<<r<<
122  "phi="<<phi<<
123  "dr="<<dr<<
124  "drphi="<<drphi<<
125  //
126  // B-Field
127  //
128  "bx="<<bx<<
129  "by="<<by<<
130  "bz="<<bz<<
131  "br="<< br<<
132  "brfi="<<brfi<<
133  // B-field integ
134  "bxi="<<bxi<<
135  "byi="<<byi<<
136  "bzi="<<bzi<<
137  "bri="<< bri<<
138  "brfii="<<brfii<<
139  "\n";
140  }
141 }
142 
143 
144 
145 Double_t AliTPCExB::GetDr(Double_t r, Double_t phi, Double_t z, Double_t bz){
148 
149  AliTPCExB *exb = Instance();
150  if (!exb) exb = AliTPCcalibDB::GetExB(bz,kFALSE);
151  if (!exb) return 0;
152  Double_t pos0[3] = {r*TMath::Cos(phi), r*TMath::Sin(phi),z};
153  Double_t pos1[3];
154  exb->Correct(pos0,pos1);
155  Double_t dx=pos1[0]-pos0[0];
156  Double_t dy=pos1[1]-pos0[1];
157  // Double_t dz=pos1[2]-pos0[2];
158  // return TMath::Sqrt(dx*dx+dy*dy);
159  Float_t dr = (dx*pos0[0]+dy*pos0[1])/r;
160  return dr;
161 }
162 
163 
164 Double_t AliTPCExB::GetDrphi(Double_t r, Double_t phi, Double_t z, Double_t bz){
166 
167  AliTPCExB *exb = Instance();
168  if (!exb) exb = AliTPCcalibDB::GetExB(bz,kFALSE);
169  if (!exb) return 0;
170  Double_t pos0[3] = {r*TMath::Cos(phi), r*TMath::Sin(phi),z};
171  Double_t pos1[3];
172  exb->Correct(pos0,pos1);
173  Double_t dphi=TMath::ATan2(pos1[1],pos1[0])-TMath::ATan2(pos0[1],pos0[0]);
174  if (dphi>TMath::Pi()) dphi-=TMath::TwoPi();
175  if (dphi<-TMath::Pi()) dphi+=TMath::TwoPi();
176  return r*dphi;
177 
178 }
179 
180 
181 Double_t AliTPCExB::GetDphi(Double_t r, Double_t phi, Double_t z, Double_t bz){
183 
184  AliTPCExB *exb = Instance();
185  if (!exb) exb = AliTPCcalibDB::GetExB(bz,kFALSE);
186  if (!exb) return 0;
187  Double_t pos0[3] = {r*TMath::Cos(phi), r*TMath::Sin(phi),z};
188  Double_t pos1[3];
189  exb->Correct(pos0,pos1);
190  Double_t dphi=TMath::ATan2(pos1[1],pos1[0])-TMath::ATan2(pos0[1],pos0[0]);
191  return dphi;
192 
193 }
194 
195 Double_t AliTPCExB::GetDz(Double_t r, Double_t phi, Double_t z, Double_t bz){
197 
198  AliTPCExB *exb = Instance();
199  if (!exb) exb = AliTPCcalibDB::GetExB(bz,kFALSE);
200  if (!exb) return 0;
201  Double_t pos0[3] = {r*TMath::Cos(phi), r*TMath::Sin(phi),z};
202  Double_t pos1[3];
203  exb->Correct(pos0,pos1);
204  Double_t dz=pos1[2]-pos0[2];
205  return dz;
206 }
207 
208 //
209 // Magnetic field
210 //
211 
212 
213 
214 
215 void AliTPCExB::RegisterField(Int_t index, AliMagF * magf){
217 
218  fgArray.AddAt(magf,index);
219 }
220 
221 
222 
223 Double_t AliTPCExB::GetBx(Double_t r, Double_t phi, Double_t z,Int_t index){
225 
226  AliMagF *mag = (AliMagF*)fgArray.At(index);
227  if (!mag) return 0;
228  Double_t xyz[3]={r*TMath::Cos(phi),r*TMath::Sin(phi),z};
229  // xyz[1]+=30;
230  Double_t bxyz[3];
231  mag->Field(xyz,bxyz);
232  return bxyz[0];
233 }
234 
235 Double_t AliTPCExB::GetBy(Double_t r, Double_t phi, Double_t z,Int_t index){
237 
238  AliMagF *mag = (AliMagF*)fgArray.At(index);
239  if (!mag) return 0;
240  Double_t xyz[3]={r*TMath::Cos(phi),r*TMath::Sin(phi),z};
241  // xyz[1]+=30;
242  Double_t bxyz[3];
243  mag->Field(xyz,bxyz);
244  return bxyz[1];
245 }
246 
247 Double_t AliTPCExB::GetBz(Double_t r, Double_t phi, Double_t z,Int_t index){
249 
250  AliMagF *mag = (AliMagF*)fgArray.At(index);
251  if (!mag) return 0;
252  Double_t xyz[3]={r*TMath::Cos(phi),r*TMath::Sin(phi),z};
253  // xyz[1]+=30;
254  Double_t bxyz[3];
255  mag->Field(xyz,bxyz);
256  return bxyz[2];
257 }
258 
259 
260 
261 Double_t AliTPCExB::GetBr(Double_t r, Double_t phi, Double_t z,Int_t index){
263 
264  AliMagF *mag = (AliMagF*)fgArray.At(index);
265  if (!mag) return 0;
266  Double_t xyz[3]={r*TMath::Cos(phi),r*TMath::Sin(phi),z};
267  //xyz[1]+=30;
268  Double_t bxyz[3];
269  mag->Field(xyz,bxyz);
270  if (r==0) return 0;
271  Double_t br = (bxyz[0]*xyz[0]+bxyz[1]*xyz[1])/r;
272  return br;
273 }
274 
275 Double_t AliTPCExB::GetBrfi(Double_t r, Double_t phi, Double_t z,Int_t index){
277 
278  AliMagF *mag = (AliMagF*)fgArray.At(index);
279  if (!mag) return 0;
280  Double_t xyz[3]={r*TMath::Cos(phi),r*TMath::Sin(phi),z};
281  //xyz[1]+=30;
282  Double_t bxyz[3];
283  mag->Field(xyz,bxyz);
284  if (r==0) return 0;
285  Double_t br = (-bxyz[0]*xyz[1]+bxyz[1]*xyz[0])/r;
286  return br;
287 }
288 
289 
290 
291 
292 Double_t AliTPCExB::GetBxI(Double_t r, Double_t phi, Double_t z,Int_t index)
293 {
294  Double_t sumf =0;
295  if (z>0 &&z<250){
296  for (Float_t zi=z;zi<250;zi+=5){
297  sumf+=GetBx(r,phi,zi,index)/GetBz(r,phi,zi,index);
298  }
299  }
300  if (z<0 &&z>-250){
301  for (Float_t zi=z;zi>-250;zi-=5){
302  sumf+=GetBx(r,phi,zi,index)/GetBz(r,phi,zi,index);
303  }
304  }
305  return sumf*5;
306 }
307 
308 Double_t AliTPCExB::GetByI(Double_t r, Double_t phi, Double_t z,Int_t index)
309 {
310  Double_t sumf =0;
311  if (z>0 &&z<250){
312  for (Float_t zi=z;zi<250;zi+=5){
313  sumf+=GetBy(r,phi,zi,index)/GetBz(r,phi,zi,index);
314  }
315  }
316  if (z<0 &&z>-250){
317  for (Float_t zi=z;zi>-250;zi-=5){
318  sumf+=GetBy(r,phi,zi,index)/GetBz(r,phi,zi,index);
319  }
320  }
321  return sumf*5;
322 }
323 
324 Double_t AliTPCExB::GetBzI(Double_t r, Double_t phi, Double_t z,Int_t index)
325 {
326  Double_t sumf =0;
327  if (z>0 &&z<250){
328  for (Float_t zi=z;zi<250;zi+=5){
329  sumf+=GetBz(r,phi,zi,index);
330  }
331  }
332  if (z<0 &&z>-250){
333  for (Float_t zi=z;zi>-250;zi-=5){
334  sumf+=GetBz(r,phi,zi,index);
335  }
336  }
337  return sumf*5;
338 }
339 
340 
341 Double_t AliTPCExB::GetBrI(Double_t r, Double_t phi, Double_t z,Int_t index)
342 {
343  Double_t sumf =0;
344  if (z>0 &&z<250){
345  for (Float_t zi=z;zi<250;zi+=5){
346  sumf+=GetBr(r,phi,zi,index)/GetBz(r,phi,zi,index);
347  }
348  }
349  if (z<0 &&z>-250){
350  for (Float_t zi=z;zi>-250;zi-=5){
351  sumf+=GetBr(r,phi,zi,index)/GetBz(r,phi,zi,index);
352  }
353  }
354  return sumf*5;
355 }
356 
357 Double_t AliTPCExB::GetBrfiI(Double_t r, Double_t phi, Double_t z,Int_t index)
358 {
359  Double_t sumf =0;
360  if (z>0 &&z<250){
361  for (Float_t zi=z;zi<250;zi+=5.){
362  sumf+=GetBrfi(r,phi,zi,index)/GetBz(r,phi,zi,index);
363  }
364  }
365  if (z<0 &&z>-250){
366  for (Float_t zi=z;zi>-250;zi-=5){
367  sumf+=GetBrfi(r,phi,zi,index)/GetBz(r,phi,zi,index);
368  }
369  }
370  return sumf*5;
371 }
372 
373 
374 Double_t AliTPCExB::Eval(Int_t type, Double_t r, Double_t phi, Double_t z){
378 
379  if (type==0) {
380  if (z>0 && fMatBrBzI0) return EvalMat(*fMatBrBzI0,r,phi,z);
381  if (z<0 && fMatBrBzI1) return EvalMat(*fMatBrBzI1,r,phi,z);
382  }
383  // brfi integral param
384  if (type==1) {
385  if (z>0 && fMatBrfiBzI0) return EvalMat(*fMatBrfiBzI0,r,phi,z);
386  if (z<0 && fMatBrfiBzI1) return EvalMat(*fMatBrfiBzI1,r,phi,z);
387  }
388  // brbz param
389  if (type==2 && fMatBrBz) return EvalMat(*fMatBrBz,r,phi,z);
390  // brfibz param
391  if (type==3 && fMatBrfiBz) return EvalMat(*fMatBrfiBz,r,phi,z);
392  return 0;
393 }
394 
395 
396 Double_t AliTPCExB::EvalMat(const TVectorD &vec, Double_t r, Double_t phi, Double_t z){
429 
430  Double_t sa = TMath::Sin(phi);
431  Double_t ca = TMath::Cos(phi);
432  Double_t sa2 = TMath::Sin(phi*2);
433  Double_t ca2 = TMath::Cos(phi*2);
434  Double_t zn = z/250.;
435  Double_t rn = r/250.;
436  Int_t ipoint=0;
437  Double_t res = vec[ipoint++];
438  res+=vec[ipoint++]*zn;
439  res+=vec[ipoint++]*rn;
440  res+=vec[ipoint++]*zn*rn;
441  res+=vec[ipoint++]*zn*zn;
442  res+=vec[ipoint++]*zn*zn*rn;
443  res+=vec[ipoint++]*zn*rn*rn;
444  //
445  res+=vec[ipoint++]*sa;
446  res+=vec[ipoint++]*ca;
447  res+=vec[ipoint++]*ca2;
448  res+=vec[ipoint++]*sa2;
449  res+=vec[ipoint++]*ca*zn;
450  res+=vec[ipoint++]*sa*zn;
451  res+=vec[ipoint++]*ca2*zn;
452  res+=vec[ipoint++]*sa2*zn;
453  res+=vec[ipoint++]*ca*zn*zn;
454  res+=vec[ipoint++]*sa*zn*zn;
455  res+=vec[ipoint++]*ca*zn*rn;
456  res+=vec[ipoint++]*sa*zn*rn;
457  return res;
458 }
459 
460 
461 
462 
463 
464 
465 
466 
467 
468 
469 
470 
471 
472 
473 /*
474 
475  AliTPCExB draw;
476  draw.RegisterField(0,new AliMagWrapCheb("Maps","Maps", 2, 1., 10., AliMagWrapCheb::k5kG));
477  draw.RegisterField(1,new AliMagFMaps("Maps","Maps", 2, 1., 10., 2));
478 
479  TF2 fbz_rz_0pi("fbz_rz_0pi","AliTPCExB::GetBz(x,0*pi,y)",0,250,-250,250);
480  fbz_rz_0pi->Draw("surf2");
481 
482  TF1 fbz_z_90_00pi("fbz_z_90_00pi","AliTPCExB::GetBz(90,0*pi,x)",-250,250);
483  TF1 fbz_z_90_05pi("fbz_z_90_05pi","AliTPCExB::GetBz(90,0.5*pi,x)",-250,250);
484  TF1 fbz_z_90_10pi("fbz_z_90_10pi","AliTPCExB::GetBz(90,1.0*pi,x)",-250,250);
485  TF1 fbz_z_90_15pi("fbz_z_90_15pi","AliTPCExB::GetBz(90,1.5*pi,x)",-250,250);
486  fbz_z_90_00pi->SetLineColor(2);
487  fbz_z_90_05pi->SetLineColor(3);
488  fbz_z_90_10pi->SetLineColor(4);
489  fbz_z_90_15pi->SetLineColor(5);
490  fbz_z_90_00pi->Draw()
491  fbz_z_90_05pi->Draw("same")
492  fbz_z_90_15pi->Draw("same")
493  fbz_z_90_10pi->Draw("same")
494 
495 
496  TF1 fbr_z_90_00pi("fbz_z_90_00pi","AliTPCExB::GetBr(90,0*pi,x)",-250,250);
497  TF1 fbr_z_90_05pi("fbz_z_90_05pi","AliTPCExB::GetBr(90,0.5*pi,x)",-250,250);
498  TF1 fbr_z_90_10pi("fbz_z_90_10pi","AliTPCExB::GetBr(90,1.0*pi,x)",-250,250);
499  TF1 fbr_z_90_15pi("fbz_z_90_15pi","AliTPCExB::GetBr(90,1.5*pi,x)",-250,250);
500  fbr_z_90_00pi->SetLineColor(2);
501  fbr_z_90_05pi->SetLineColor(3);
502  fbr_z_90_10pi->SetLineColor(4);
503  fbr_z_90_15pi->SetLineColor(5);
504  fbr_z_90_00pi->Draw()
505  fbr_z_90_05pi->Draw("same")
506  fbr_z_90_15pi->Draw("same")
507  fbr_z_90_10pi->Draw("same")
508 
509  //
510  TF2 fbz_xy_0z("fbz_xy_0z","AliTPCExB::GetBz(sqrt(x^2+y^2),atan2(y,x),0)",-250,250,-250,250);
511  fbz_xy_0z.SetNpy(100);
512  fbz_xy_0z.SetNpx(100);
513  fbz_xy_0z->Draw("colz");
514  //
515  TF2 fbz_xy_250z("fbz_xy_250z","AliTPCExB::GetBz(sqrt(x^2+y^2),atan2(y,x),250)",-250,250,-250,250);
516  fbz_xy_250z.SetNpy(100);
517  fbz_xy_250z.SetNpx(100)
518  fbz_xy_250z->Draw("colz");
519  //
520  TF2 fbz_xy_m250z("fbz_xy_m250z","AliTPCExB::GetBz(sqrt(x^2+y^2),atan2(y,x),-250)",-250,250,-250,250);
521  fbz_xy_m250z.SetNpy(100);
522  fbz_xy_m250z.SetNpx(100)
523  fbz_xy_m250z->Draw("colz");
524  //
525 
526 
527 
528 
529 
530 */
531 
532 
static Double_t GetDz(Double_t r, Double_t phi, Double_t z, Double_t bz=5)
Definition: AliTPCExB.cxx:195
TVectorD * fMatBrBz
param matrix Br/Bz
Definition: AliTPCExB.h:54
TVectorD * fMatBrBzI1
param matrix Br/Bz integral z<0
Definition: AliTPCExB.h:57
#define TObjArray
TVectorD vec
Definition: AnalyzeLaser.C:8
static Double_t GetDrphi(Double_t r, Double_t phi, Double_t z, Double_t bz=5)
Definition: AliTPCExB.cxx:164
static Double_t GetBxI(Double_t r, Double_t phi, Double_t z, Int_t index=0)
Definition: AliTPCExB.cxx:292
static Double_t GetByI(Double_t r, Double_t phi, Double_t z, Int_t index=0)
Definition: AliTPCExB.cxx:308
virtual void Field(const Double_t *x, Double_t *b)
Definition: AliMagF.cxx:280
static Double_t GetBrI(Double_t r, Double_t phi, Double_t z, Int_t index=0)
Definition: AliTPCExB.cxx:341
AliTPCExB & operator=(const AliTPCExB &exb)
Definition: AliTPCExB.cxx:65
AliTPCExB * GetExB() const
Definition: AliTPCcalibDB.h:64
static Double_t GetBrfi(Double_t r, Double_t phi, Double_t z, Int_t index=0)
Definition: AliTPCExB.cxx:275
static Double_t GetBy(Double_t r, Double_t phi, Double_t z, Int_t index=0)
Definition: AliTPCExB.cxx:235
TString fileName(const char *dir, int runNumber, const char *da, int i, const char *type)
static TObjArray fgArray
! array of magnetic fields
Definition: AliTPCExB.h:62
static AliTPCExB * fgInstance
! Instance of this class (singleton implementation)
Definition: AliTPCExB.h:61
void TestExB(const char *fileName)
Definition: AliTPCExB.cxx:75
static Double_t GetBz(Double_t r, Double_t phi, Double_t z, Int_t index=0)
Definition: AliTPCExB.cxx:247
static Double_t GetBr(Double_t r, Double_t phi, Double_t z, Int_t index=0)
Definition: AliTPCExB.cxx:261
static Double_t GetBrfiI(Double_t r, Double_t phi, Double_t z, Int_t index=0)
Definition: AliTPCExB.cxx:357
static void RegisterField(Int_t index, AliMagF *magf)
Definition: AliTPCExB.cxx:215
Double_t Eval(Int_t type, Double_t r, Double_t phi, Double_t z)
Definition: AliTPCExB.cxx:374
static Double_t GetDr(Double_t r, Double_t phi, Double_t z, Double_t bz=5)
Definition: AliTPCExB.cxx:145
static Double_t GetDphi(Double_t r, Double_t phi, Double_t z, Double_t bz=5)
Definition: AliTPCExB.cxx:181
TVectorD * fMatBrfiBzI1
param matrix Br/Bz integral z<0
Definition: AliTPCExB.h:59
virtual void Correct(const Double_t *position, Double_t *corrected)=0
TVectorD * fMatBrfiBz
param matrix Br/Bz
Definition: AliTPCExB.h:55
static AliMagF::BMap_t mag
class TVectorT< Double_t > TVectorD
TVectorD * fMatBrfiBzI0
param matrix Br/Bz integral z>0
Definition: AliTPCExB.h:58
void res(Char_t i)
Definition: Resolution.C:2
TVectorD * fMatBrBzI0
param matrix Br/Bz integral z>0
Definition: AliTPCExB.h:56
Abstract class for ExB effect parameterization.
Definition: AliTPCExB.h:10
static Double_t EvalMat(const TVectorD &vec, Double_t r, Double_t phi, Double_t z)
Definition: AliTPCExB.cxx:396
static Double_t GetBx(Double_t r, Double_t phi, Double_t z, Int_t index=0)
Definition: AliTPCExB.cxx:223
static AliTPCExB * Instance()
Definition: AliTPCExB.h:30
static Double_t GetBzI(Double_t r, Double_t phi, Double_t z, Int_t index=0)
Definition: AliTPCExB.cxx:324