]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCExB.cxx
Updated flags for low flux case (A. Dainese)
[u/mrichter/AliRoot.git] / TPC / AliTPCExB.cxx
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
9 //
10 // Abstract class for ExB effect parameterization
11 // 
12 //
13 // 
14 // The ExB correction map is stored in the calib DB
15 // The lookup can be dumped to the tree:
16 /*
17
18    //
19   char *storage = "local://OCDBres"
20   Int_t RunNumber=0;
21   AliCDBManager::Instance()->SetDefaultStorage(storage);
22   AliCDBManager::Instance()->SetRun(RunNumber) 
23   AliTPCExBFirst * exb = AliTPCcalibDB::Instance()->GetExB();
24   //
25   // See example macro $ALICE_ROOT/TPC/macros/AliTPCExBdraw.C 
26   //
27   .L $ALICE_ROOT/TPC/macros/AliTPCExBdraw.C 
28   Draw(0)
29
30
31
32
33 */
34
35 AliTPCExB* AliTPCExB::fgInstance = 0;
36
37 TObjArray   AliTPCExB::fgArray;
38
39
40
41 ClassImp(AliTPCExB)
42
43 AliTPCExB::AliTPCExB():
44   TObject(),
45   fMatBrBz(0),       //param matrix Br/Bz
46   fMatBrfiBz(0),     //param matrix Br/Bz
47   fMatBrBzI0(0),     //param matrix Br/Bz integral  z>0 
48   fMatBrBzI1(0),     //param matrix Br/Bz integral  z<0 
49   fMatBrfiBzI0(0),   //param matrix Br/Bz integral  z>0 
50   fMatBrfiBzI1(0)    //param matrix Br/Bz integral  z<0
51 {
52   //
53   // default constructor
54   //
55 }
56
57 AliTPCExB::AliTPCExB(const AliTPCExB& exb):
58   TObject(exb),
59   fMatBrBz(new TVectorD(*(exb.fMatBrBz))),       //param matrix Br/Bz
60   fMatBrfiBz(new TVectorD(*(exb.fMatBrfiBz))),     //param matrix Br/Bz
61   fMatBrBzI0(new TVectorD(*(exb.fMatBrBzI0))),     //param matrix Br/Bz integral  z>0 
62   fMatBrBzI1(new TVectorD(*(exb.fMatBrBzI1))),     //param matrix Br/Bz integral  z<0 
63   fMatBrfiBzI0(new TVectorD(*(exb.fMatBrfiBzI0))),   //param matrix Br/Bz integral  z>0 
64   fMatBrfiBzI1(new TVectorD(*(exb.fMatBrfiBzI1)))    //param matrix Br/Bz integral  z<0
65 {
66   //
67   // copy constructor
68   //
69 }
70
71 AliTPCExB& AliTPCExB::operator=(const AliTPCExB &/*exb*/)
72 {
73   //
74   // Dummy  assignment
75   //
76   return *this;
77 }
78
79
80
81
82 void AliTPCExB::TestExB(const char* fileName) {
83   //
84   // Test ExB  - 
85   // Dump the filed and corrections to the tree in file fileName  
86   //
87   // 
88   TTreeSRedirector ts(fileName);
89   Double_t x[3];
90   for (x[0]=-250.;x[0]<=250.;x[0]+=10.)
91     for (x[1]=-250.;x[1]<=250.;x[1]+=10.)
92       for (x[2]=-250.;x[2]<=250.;x[2]+=20.) {
93         Double_t r=TMath::Sqrt(x[0]*x[0]+x[1]*x[1]);
94         if (r<20) continue;
95         if (r>260) continue;
96         Double_t z = x[2];
97         Double_t d[3];
98         Correct(x,d);
99         Double_t rd=TMath::Sqrt(d[0]*d[0]+d[1]*d[1]);
100         Double_t dr=r-rd;
101         Double_t phi=TMath::ATan2(x[1],x[0]);
102         Double_t phid=TMath::ATan2(d[1],d[0]);
103         Double_t dphi=phi-phid;
104         if (dphi<0.) dphi+=TMath::TwoPi();
105         if (dphi>TMath::Pi()) dphi=TMath::TwoPi()-dphi;
106         Double_t drphi=r*dphi;
107         Double_t dx=x[0]-d[0];
108         Double_t dy=x[1]-d[1];
109         Double_t dz=x[2]-d[2];
110         //
111         Double_t bx    = GetBx(r,phi,z,0);
112         Double_t by    = GetBy(r,phi,z,0);
113         Double_t bz    = GetBz(r,phi,z,0);
114         Double_t br    = GetBr(r,phi,z,0);
115         Double_t brfi  = GetBrfi(r,phi,z,0);
116         //
117         Double_t bxi    = GetBxI(r,phi,z,0);
118         Double_t byi    = GetByI(r,phi,z,0);
119         Double_t bzi    = GetBzI(r,phi,z,0);
120         Double_t bri    = GetBrI(r,phi,z,0);
121         Double_t brfii  = GetBrfiI(r,phi,z,0);
122
123         ts<<"positions"<<
124           "x0="<<x[0]<<
125           "x1="<<x[1]<<
126           "x2="<<x[2]<<
127           "dx="<<dx<<
128           "dy="<<dy<<
129           "dz="<<dz<<
130           "r="<<r<<
131           "phi="<<phi<<
132           "dr="<<dr<<
133           "drphi="<<drphi<<
134           //
135           // B-Field
136           //
137           "bx="<<bx<<
138           "by="<<by<<
139           "bz="<<bz<<
140           "br="<< br<<
141           "brfi="<<brfi<<
142           // B-field integ
143           "bxi="<<bxi<<
144           "byi="<<byi<<
145           "bzi="<<bzi<<
146           "bri="<< bri<<
147           "brfii="<<brfii<<
148           "\n";
149       }
150 }
151
152
153
154 Double_t AliTPCExB::GetDr(Double_t r, Double_t phi, Double_t z, Double_t bz){
155   //
156   // Static function
157   // Posibble to us it for visualization 
158   // 
159   //
160   AliTPCExB *exb = AliTPCcalibDB::GetExB(bz,kFALSE);
161   if (!exb) return 0;
162   Double_t pos0[3] = {r*TMath::Cos(phi), r*TMath::Sin(phi),z};
163   Double_t pos1[3];
164   exb->Correct(pos0,pos1);
165   Double_t dx=pos1[0]-pos0[0];
166   Double_t dy=pos1[1]-pos0[1];
167   //  Double_t dz=pos1[2]-pos0[2];
168   // return TMath::Sqrt(dx*dx+dy*dy);  
169   Float_t dr = (dx*pos0[0]+dy*pos0[1])/r;
170   return dr;
171 }
172
173
174 Double_t AliTPCExB::GetDrphi(Double_t r, Double_t phi, Double_t z, Double_t bz){
175   //
176   //
177   //
178   AliTPCExB *exb = AliTPCcalibDB::GetExB(bz,kFALSE);
179   if (!exb) return 0;
180   Double_t pos0[3] = {r*TMath::Cos(phi), r*TMath::Sin(phi),z};
181   Double_t pos1[3];
182   exb->Correct(pos0,pos1);
183   Double_t dphi=TMath::ATan2(pos1[1],pos1[0])-TMath::ATan2(pos0[1],pos0[0]);
184   return r*dphi;  
185
186 }
187
188
189 Double_t AliTPCExB::GetDphi(Double_t r, Double_t phi, Double_t z, Double_t bz){
190   //
191   //
192   // 
193   AliTPCExB *exb = AliTPCcalibDB::GetExB(bz,kFALSE);
194   if (!exb) return 0;
195   Double_t pos0[3] = {r*TMath::Cos(phi), r*TMath::Sin(phi),z};
196   Double_t pos1[3];
197   exb->Correct(pos0,pos1);
198   Double_t dphi=TMath::ATan2(pos1[1],pos1[0])-TMath::ATan2(pos0[1],pos0[0]);
199   return dphi;  
200
201 }
202
203 Double_t AliTPCExB::GetDz(Double_t r, Double_t phi, Double_t z, Double_t bz){
204   //
205   //
206   //
207   AliTPCExB *exb = AliTPCcalibDB::GetExB(bz,kFALSE);
208   if (!exb) return 0;
209   Double_t pos0[3] = {r*TMath::Cos(phi), r*TMath::Sin(phi),z};
210   Double_t pos1[3];
211   exb->Correct(pos0,pos1);
212   Double_t dz=pos1[2]-pos0[2];
213   return dz;  
214 }
215
216 //
217 // Magnetic field
218 //
219
220
221
222
223 void AliTPCExB::RegisterField(Int_t index, AliMagF * magf){
224   //
225   // add the filed to the list
226   //
227   fgArray.AddAt(magf,index);
228 }
229
230
231
232 Double_t AliTPCExB::GetBx(Double_t r, Double_t phi, Double_t z,Int_t index){
233   //
234   // 
235   //
236   AliMagF *mag = (AliMagF*)fgArray.At(index);
237   if (!mag) return 0;
238   Float_t xyz[3]={r*TMath::Cos(phi),r*TMath::Sin(phi),z};
239   //  xyz[1]+=30;
240   Float_t bxyz[3];
241   mag->Field(xyz,bxyz);
242   return bxyz[0];
243 }  
244
245 Double_t AliTPCExB::GetBy(Double_t r, Double_t phi, Double_t z,Int_t index){
246   //
247   // 
248   //
249   AliMagF *mag = (AliMagF*)fgArray.At(index);
250   if (!mag) return 0;
251   Float_t xyz[3]={r*TMath::Cos(phi),r*TMath::Sin(phi),z};
252   //  xyz[1]+=30;
253   Float_t bxyz[3];
254   mag->Field(xyz,bxyz);
255   return bxyz[1];
256 }  
257
258 Double_t AliTPCExB::GetBz(Double_t r, Double_t phi, Double_t z,Int_t index){
259   //
260   // 
261   //
262   AliMagF *mag = (AliMagF*)fgArray.At(index);
263   if (!mag) return 0;
264   Float_t xyz[3]={r*TMath::Cos(phi),r*TMath::Sin(phi),z};
265   //  xyz[1]+=30;
266   Float_t bxyz[3];
267   mag->Field(xyz,bxyz);
268   return bxyz[2];
269 }  
270
271
272
273 Double_t AliTPCExB::GetBr(Double_t r, Double_t phi, Double_t z,Int_t index){
274   //
275   // 
276   //
277   AliMagF *mag = (AliMagF*)fgArray.At(index);
278   if (!mag) return 0;
279   Float_t xyz[3]={r*TMath::Cos(phi),r*TMath::Sin(phi),z};
280   //xyz[1]+=30;
281   Float_t bxyz[3];
282   mag->Field(xyz,bxyz);
283   if (r==0) return 0;
284   Float_t br = (bxyz[0]*xyz[0]+bxyz[1]*xyz[1])/r;
285   return br;
286 }  
287
288 Double_t AliTPCExB::GetBrfi(Double_t r, Double_t phi, Double_t z,Int_t index){
289   //
290   // 
291   //
292   AliMagF *mag = (AliMagF*)fgArray.At(index);
293   if (!mag) return 0;
294   Float_t xyz[3]={r*TMath::Cos(phi),r*TMath::Sin(phi),z};
295   //xyz[1]+=30;
296   Float_t bxyz[3];
297   mag->Field(xyz,bxyz);
298   if (r==0) return 0;
299   Float_t br = (-bxyz[0]*xyz[1]+bxyz[1]*xyz[0])/r;
300   return br;
301 }  
302
303
304
305
306 Double_t AliTPCExB::GetBxI(Double_t r, Double_t phi, Double_t z,Int_t index)
307 {
308   Double_t sumf  =0;
309   if (z>0 &&z<250){
310     for (Float_t zi=z;zi<250;zi+=5){
311       sumf+=GetBx(r,phi,zi,index)/GetBz(r,phi,zi,index);
312     }
313   }
314   if (z<0 &&z>-250){
315     for (Float_t zi=z;zi>-250;zi-=5){
316       sumf+=GetBx(r,phi,zi,index)/GetBz(r,phi,zi,index);
317     }
318   }
319   return sumf*5;
320 }
321
322 Double_t AliTPCExB::GetByI(Double_t r, Double_t phi, Double_t z,Int_t index)
323 {
324   Double_t sumf =0;
325   if (z>0 &&z<250){
326     for (Float_t zi=z;zi<250;zi+=5){
327       sumf+=GetBy(r,phi,zi,index)/GetBz(r,phi,zi,index);
328     }
329   }
330   if (z<0 &&z>-250){
331     for (Float_t zi=z;zi>-250;zi-=5){
332       sumf+=GetBy(r,phi,zi,index)/GetBz(r,phi,zi,index);
333     }
334   }
335   return sumf*5;
336 }
337
338 Double_t AliTPCExB::GetBzI(Double_t r, Double_t phi, Double_t z,Int_t index)
339 {
340   Double_t sumf =0;
341   if (z>0 &&z<250){
342     for (Float_t zi=z;zi<250;zi+=5){
343       sumf+=GetBz(r,phi,zi,index);
344     }
345   }
346   if (z<0 &&z>-250){
347     for (Float_t zi=z;zi>-250;zi-=5){
348       sumf+=GetBz(r,phi,zi,index);
349     }
350   }
351   return sumf*5;
352 }
353
354
355 Double_t AliTPCExB::GetBrI(Double_t r, Double_t phi, Double_t z,Int_t index)
356 {
357   Double_t sumf =0;
358   if (z>0 &&z<250){
359     for (Float_t zi=z;zi<250;zi+=5){
360       sumf+=GetBr(r,phi,zi,index)/GetBz(r,phi,zi,index);
361     }
362   }
363   if (z<0 &&z>-250){
364     for (Float_t zi=z;zi>-250;zi-=5){
365       sumf+=GetBr(r,phi,zi,index)/GetBz(r,phi,zi,index);
366     }
367   }
368   return sumf*5;
369 }
370
371 Double_t AliTPCExB::GetBrfiI(Double_t r, Double_t phi, Double_t z,Int_t index)
372 {
373   Double_t sumf =0;
374   if (z>0 &&z<250){
375     for (Float_t zi=z;zi<250;zi+=5.){
376       sumf+=GetBrfi(r,phi,zi,index)/GetBz(r,phi,zi,index);
377     }
378   }
379   if (z<0 &&z>-250){
380     for (Float_t zi=z;zi>-250;zi-=5){
381       sumf+=GetBrfi(r,phi,zi,index)/GetBz(r,phi,zi,index);
382     }
383   }
384   return sumf*5;
385 }
386
387
388 Double_t AliTPCExB::Eval(Int_t type, Double_t r, Double_t phi, Double_t z){
389   //
390   // Evaluate parameterization
391   //
392   // br integral param 
393   if (type==0) {
394     if (z>0 && fMatBrBzI0) return EvalMat(*fMatBrBzI0,r,phi,z);
395     if (z<0 && fMatBrBzI1) return EvalMat(*fMatBrBzI1,r,phi,z);
396   }
397   // brfi integral param   
398   if (type==1) {
399     if (z>0 && fMatBrfiBzI0) return EvalMat(*fMatBrfiBzI0,r,phi,z);
400     if (z<0 && fMatBrfiBzI1) return EvalMat(*fMatBrfiBzI1,r,phi,z);
401   }
402   // brbz param
403   if (type==2 && fMatBrBz) return EvalMat(*fMatBrBz,r,phi,z);
404   // brfibz param
405   if (type==3 && fMatBrfiBz) return EvalMat(*fMatBrfiBz,r,phi,z);
406   return 0;
407 }
408
409
410 Double_t AliTPCExB::EvalMat(TVectorD &vec, Double_t r, Double_t phi, Double_t z){
411   //
412   // Evaluate taylor expansion in r,phi,z
413   //
414   // Variables  
415   //tree->SetAlias("sa","sin(phi+0.0)");
416   //tree->SetAlias("ca","cos(phi+0.0)");
417   //tree->SetAlias("sa2","sin(phi*2+0.0)");
418   //tree->SetAlias("ca2","cos(phi*2+0.0)");
419   //tree->SetAlias("zn","(x2/250.)");
420   //tree->SetAlias("rn","(r/250.)")
421   //   TString fstringSym="";
422   //   //  
423   //   fstringSym+="zn++";
424   //   fstringSym+="rn++";
425   //   fstringSym+="zn*rn++";
426   //   fstringSym+="zn*zn++";
427   //   fstringSym+="zn*zn*rn++";
428   //   fstringSym+="zn*rn*rn++";
429   //   //
430   //   fstringSym+="sa++";
431   //   fstringSym+="ca++";  
432   //   fstringSym+="ca2++";
433   //   fstringSym+="sa2++";
434   //   fstringSym+="ca*zn++";
435   //   fstringSym+="sa*zn++";
436   //   fstringSym+="ca2*zn++";
437   //   fstringSym+="sa2*zn++";
438   //   fstringSym+="ca*zn*zn++";
439   //   fstringSym+="sa*zn*zn++";
440   //   fstringSym+="ca*zn*rn++";
441   //   fstringSym+="sa*zn*rn++";
442
443
444   Double_t sa  = TMath::Sin(phi);
445   Double_t ca  = TMath::Cos(phi);
446   Double_t sa2 = TMath::Sin(phi*2);
447   Double_t ca2 = TMath::Cos(phi*2);
448   Double_t zn  = z/250.;
449   Double_t rn  = r/250.;
450   Int_t  ipoint=0;
451   Double_t res = vec[ipoint++];
452   res+=vec[ipoint++]*zn;
453   res+=vec[ipoint++]*rn;
454   res+=vec[ipoint++]*zn*rn;
455   res+=vec[ipoint++]*zn*zn;
456   res+=vec[ipoint++]*zn*zn*rn;
457   res+=vec[ipoint++]*zn*rn*rn;
458   //
459   res+=vec[ipoint++]*sa;
460   res+=vec[ipoint++]*ca;  
461   res+=vec[ipoint++]*ca2;
462   res+=vec[ipoint++]*sa2;
463   res+=vec[ipoint++]*ca*zn;
464   res+=vec[ipoint++]*sa*zn;
465   res+=vec[ipoint++]*ca2*zn;
466   res+=vec[ipoint++]*sa2*zn;
467   res+=vec[ipoint++]*ca*zn*zn;
468   res+=vec[ipoint++]*sa*zn*zn;
469   res+=vec[ipoint++]*ca*zn*rn;
470   res+=vec[ipoint++]*sa*zn*rn;
471   return res;
472 }
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487 /*
488   
489  AliTPCExB draw;
490  draw.RegisterField(0,new AliMagWrapCheb("Maps","Maps", 2, 1., 10., AliMagWrapCheb::k5kG));
491  draw.RegisterField(1,new AliMagFMaps("Maps","Maps", 2, 1., 10., 2));
492
493  TF2 fbz_rz_0pi("fbz_rz_0pi","AliTPCExB::GetBz(x,0*pi,y)",0,250,-250,250);
494  fbz_rz_0pi->Draw("surf2");
495  
496  TF1 fbz_z_90_00pi("fbz_z_90_00pi","AliTPCExB::GetBz(90,0*pi,x)",-250,250);
497   TF1 fbz_z_90_05pi("fbz_z_90_05pi","AliTPCExB::GetBz(90,0.5*pi,x)",-250,250);
498   TF1 fbz_z_90_10pi("fbz_z_90_10pi","AliTPCExB::GetBz(90,1.0*pi,x)",-250,250);
499   TF1 fbz_z_90_15pi("fbz_z_90_15pi","AliTPCExB::GetBz(90,1.5*pi,x)",-250,250);
500   fbz_z_90_00pi->SetLineColor(2);
501   fbz_z_90_05pi->SetLineColor(3);
502   fbz_z_90_10pi->SetLineColor(4);
503   fbz_z_90_15pi->SetLineColor(5);
504   fbz_z_90_00pi->Draw()
505   fbz_z_90_05pi->Draw("same")
506   fbz_z_90_15pi->Draw("same")
507   fbz_z_90_10pi->Draw("same")
508   
509
510   TF1 fbr_z_90_00pi("fbz_z_90_00pi","AliTPCExB::GetBr(90,0*pi,x)",-250,250);
511   TF1 fbr_z_90_05pi("fbz_z_90_05pi","AliTPCExB::GetBr(90,0.5*pi,x)",-250,250);
512   TF1 fbr_z_90_10pi("fbz_z_90_10pi","AliTPCExB::GetBr(90,1.0*pi,x)",-250,250);
513   TF1 fbr_z_90_15pi("fbz_z_90_15pi","AliTPCExB::GetBr(90,1.5*pi,x)",-250,250);
514   fbr_z_90_00pi->SetLineColor(2);
515   fbr_z_90_05pi->SetLineColor(3);
516   fbr_z_90_10pi->SetLineColor(4);
517   fbr_z_90_15pi->SetLineColor(5);
518   fbr_z_90_00pi->Draw()
519   fbr_z_90_05pi->Draw("same")
520   fbr_z_90_15pi->Draw("same")
521   fbr_z_90_10pi->Draw("same")
522
523   //
524   TF2 fbz_xy_0z("fbz_xy_0z","AliTPCExB::GetBz(sqrt(x^2+y^2),atan2(y,x),0)",-250,250,-250,250);
525   fbz_xy_0z.SetNpy(100);
526   fbz_xy_0z.SetNpx(100);
527   fbz_xy_0z->Draw("colz");
528   //
529   TF2 fbz_xy_250z("fbz_xy_250z","AliTPCExB::GetBz(sqrt(x^2+y^2),atan2(y,x),250)",-250,250,-250,250);
530   fbz_xy_250z.SetNpy(100);
531   fbz_xy_250z.SetNpx(100)
532   fbz_xy_250z->Draw("colz");
533   //
534    TF2 fbz_xy_m250z("fbz_xy_m250z","AliTPCExB::GetBz(sqrt(x^2+y^2),atan2(y,x),-250)",-250,250,-250,250);
535   fbz_xy_m250z.SetNpy(100);
536   fbz_xy_m250z.SetNpx(100)
537   fbz_xy_m250z->Draw("colz");
538   //
539
540
541
542
543
544 */
545
546