]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HBTAN/AliHBTTrackPoints.cxx
f7ec1a042dc0afda38ad323ea8422e58a6ddd24b
[u/mrichter/AliRoot.git] / HBTAN / AliHBTTrackPoints.cxx
1 #include "AliHBTTrackPoints.h"
2 //_________________________________
3 ////////////////////////////////////////////////////////////
4 //                                                        //
5 // class AliHBTTrackPoints                                //
6 //                                                        //
7 // used by Anti-Merging cut                               //
8 // contains set of poits the lay on track trajectory      //
9 // according to reconstructed track parameters -          //
10 // NOT CLUSTERS POSITIONS!!!                              //
11 // Anti-Merging cut is applied only on tracks coming from //
12 // different events (that are use to fill deniminators)   //
13 //                                                        //
14 ////////////////////////////////////////////////////////////
15 #include "AliTrackReference.h"
16
17 #include "AliTPCtrack.h"
18 #include "AliESDtrack.h"
19
20 #include <TMath.h>
21 #include <TClonesArray.h>
22
23 ClassImp(AliHBTTrackPoints)
24
25 Int_t AliHBTTrackPoints::fgDebug = 0;
26 AliHBTTrackPoints::AliHBTTrackPoints():
27  fN(0),
28  fX(0x0),
29  fY(0x0),
30  fZ(0x0)
31 {
32   //constructor
33 }
34 /***************************************************************/
35
36 AliHBTTrackPoints::AliHBTTrackPoints(Int_t n, AliESDtrack* track, Float_t mf, Float_t dr, Float_t r0):
37  fN(n),
38  fX(new Float_t[fN]),
39  fY(new Float_t[fN]),
40  fZ(new Float_t[fN])
41 {
42   //constructor
43   //mf - magnetic field in kG - needed to calculated curvature out of Pt
44   //r0 - starting radius
45   //dr - calculate points every dr cm, default every 30cm
46   if (track == 0x0)
47    {
48      Error("AliHBTTrackPoints","ESD track is null");
49      fN = 0;
50      delete [] fX;
51      delete [] fY;
52      delete [] fZ;
53      fX = fY = fZ = 0x0;
54      return;
55    }
56
57   Double_t x;
58   Double_t par[5];
59   track->GetInnerExternalParameters(x,par);     //get properties of the track
60   if (x == 0)
61    {
62      Error("AliHBTTrackPoints","This ESD track does not contain TPC information");
63      return;
64    }
65
66   if (mf == 0.0)
67    {
68      Error("AliHBTTrackPoints","Zero Magnetic field passed as parameter.");
69      return;
70    }
71    
72   Double_t alpha = track->GetInnerAlpha();
73   Double_t cc = 1000./0.299792458/mf;//conversion constant
74   Double_t c=par[4]/cc;
75   
76   MakePoints(dr,r0,x,par,c,alpha);
77   
78 }
79 /***************************************************************/
80
81 AliHBTTrackPoints::AliHBTTrackPoints(Int_t n, AliTPCtrack* track, Float_t dr, Float_t r0):
82  fN(n),
83  fX(new Float_t[fN]),
84  fY(new Float_t[fN]),
85  fZ(new Float_t[fN])
86 {
87   //constructor
88   //r0 starting radius
89   //dr - calculate points every dr cm, default every 30cm
90   if (track == 0x0)
91    {
92      Error("AliHBTTrackPoints","TPC track is null");
93      fN = 0;
94      delete [] fX;
95      delete [] fY;
96      delete [] fZ;
97      fX = fY = fZ = 0x0;
98      return;
99    }
100   track->PropagateTo(r0);
101
102  //* This formation is now fixed in the following way:                         *
103  //*      external param0:   local Y-coordinate of a track (cm)                *
104  //*      external param1:   local Z-coordinate of a track (cm)                *
105  //*      external param2:   local sine of the track momentum azimuth angle    *
106  //*      external param3:   tangent of the track momentum dip angle           *
107  //*      external param4:   1/pt (1/(GeV/c))                                  *
108     
109   Double_t x = 0;
110   Double_t par[5];
111   track->GetExternalParameters(x,par);     //get properties of the track
112   
113   Double_t alpha = track->GetAlpha();
114   Double_t c=track->GetC();
115   MakePoints(dr,r0,x,par,c,alpha);
116 }  
117
118 void AliHBTTrackPoints::MakePoints( Float_t dr, Float_t r0, Double_t x, Double_t* par, Double_t c, Double_t alpha)
119 {
120   //Calculates points starting at radius r0
121   //spacing every dr (in radial direction)
122   // according to track parameters
123   // x - position in sector local reference frame. x i parallel to R and sector is symmetric with respect to x
124   // par - track external parameters; array with 5 elements;  look at AliTPCtrack.h or AliESDtrack.h for their meaning
125   // c - track curvature
126   // alpha - sector's rotation angle (phi) == angle needed for local to global transformation
127   
128   Double_t y = par[0];
129   Double_t z0 = par[1];
130       
131   Double_t phi0local = TMath::ATan2(y,x);
132   Double_t phi0global = phi0local + alpha;
133   
134   if (phi0local<0) phi0local+=2*TMath::Pi();
135   if (phi0local>=2.*TMath::Pi()) phi0local-=2*TMath::Pi();
136   
137   if (phi0global<0) phi0global+=2*TMath::Pi();
138   if (phi0global>=2.*TMath::Pi()) phi0global-=2*TMath::Pi();
139   
140   Double_t r = TMath::Hypot(x,y);
141   
142
143   if (GetDebug() > 9) 
144     Info("AliHBTTrackPoints","Radius0 %f, Real Radius %f",r0,r); 
145   
146   if (GetDebug() > 5) 
147     Info("AliHBTTrackPoints","Phi Global at first padraw %f, Phi locat %f",phi0global,phi0local);
148     
149   Double_t eta = x*c - par[2] ;//par[2] = fX*C - eta; eta==fP2 ; C==fP4
150   
151   //this calculattions are assuming new (current) model 
152   Double_t tmp = par[2];
153   tmp = 1. - tmp*tmp;
154   tmp = c*y + TMath::Sqrt(tmp);
155   Double_t dca=(TMath::Hypot(eta,tmp) - 1. )/TMath::Abs(c);
156
157   //Here is old model Cold=Cnew/2.
158   Double_t dcasq = dca*dca;
159   Double_t c2 = c/2.;
160   Double_t cst1 = (1.+c2*dca)*dca;//first constant
161   Double_t cst2 = 1. + 2.*c2*dca;//second constant
162       
163   Double_t factorPhi0 = TMath::ASin((c2*r + cst1/r)/cst2);
164   Double_t factorZ0 = TMath::ASin(c2*TMath::Sqrt((r*r-dcasq)/cst2))*par[3]/c2;
165  
166   for(Int_t i = 0; i<fN; i++)
167    {
168      Double_t rc = r0 + i*dr;
169      Double_t factorPhi = TMath::ASin( (c2*rc + cst1/rc)/cst2 );//factor phi od rc
170      Double_t phi = phi0global + factorPhi - factorPhi0;
171      
172      Double_t factorZ = TMath::ASin(c2*TMath::Sqrt((rc*rc-dcasq)/cst2))*par[3]/c2;
173      fZ[i] = z0 + factorZ - factorZ0;
174      fX[i] = rc*TMath::Cos(phi);
175      fY[i] = rc*TMath::Sin(phi);
176      
177      if ( GetDebug() > 2 )
178       {
179         Info("AliHBTTrackPoints","X %f Y %f Z %f R asked %f R obtained %f",
180              fX[i],fY[i],fZ[i],rc,TMath::Hypot(fX[i],fY[i]));
181       }
182    }
183 }
184 /***************************************************************/
185
186 AliHBTTrackPoints::~AliHBTTrackPoints()
187 {
188   //destructor
189   delete [] fX;
190   delete [] fY;
191   delete [] fZ;
192 }
193 /***************************************************************/
194
195 void AliHBTTrackPoints::PositionAt(Int_t n, Float_t &x,Float_t &y,Float_t &z)
196 {
197   //returns position at point n
198   if ((n<0) || (n>fN))
199    {
200      Error("PositionAt","Point %d out of range",n);
201      return;
202    }
203   
204   x = fX[n];
205   y = fY[n];
206   z = fZ[n];
207   if ( GetDebug() > 1 )
208     {
209       Info("AliHBTTrackPoints","n %d; X %f; Y %f; Z %f",n,x,y,z);
210     }
211 }
212 /***************************************************************/
213
214 Double_t AliHBTTrackPoints::AvarageDistance(const AliHBTTrackPoints& tr)
215 {
216   //returns the aritmethic avarage distance between two tracks
217   if (fN != tr.fN)
218    {
219      Error("AvarageDistance","Number of points is not equal");
220      return -1;
221    }
222   if ( (fN <= 0) || (tr.fN <=0) )
223    {
224      Error("AvarageDistance","One of tracks is empty");
225      return -1;
226    }
227    
228   Double_t sum = 0;
229   for (Int_t i = 0; i<fN; i++)
230    {
231      if (GetDebug()>9)
232       {
233 //       Float_t r1sq = fX[i]*fX[i]+fY[i]*fY[i];
234 //       Float_t r2sq = tr.fX[i]*tr.fX[i]+tr.fY[i]*tr.fY[i];
235        Float_t r1sq = TMath::Hypot(fX[i],fY[i]);
236        Float_t r2sq = TMath::Hypot(tr.fX[i],tr.fY[i]);
237        Info("AvarageDistance","radii: %f %f",r1sq,r2sq);
238       } 
239
240       
241      Double_t dx = fX[i]-tr.fX[i];
242      Double_t dy = fY[i]-tr.fY[i];
243      Double_t dz = fZ[i]-tr.fZ[i];
244      sum+=TMath::Sqrt(dx*dx + dy*dy + dz*dz);
245      
246      if (GetDebug()>0)
247       {
248        Info("AvarageDistance","Diff: x ,y z: %f , %f, %f",dx,dy,dz);
249        Info("AvarageDistance","xxyyzz %f %f %f %f %f %f",
250             fX[i],tr.fX[i],fY[i],tr.fY[i],fZ[i],tr.fZ[i]);
251       } 
252    }
253    
254   Double_t retval = sum/((Double_t)fN);
255   if ( GetDebug() )
256     {
257       Info("AvarageDistance","Avarage distance is %f.",retval);
258     }
259   return retval;
260 }
261 /***************************************************************/
262 /***************************************************************/
263 /***************************************************************/
264 /***************************************************************/
265 /***************************************************************/
266 /***************************************************************/
267
268 #include "AliRun.h"
269 #include "AliESD.h"
270 #include "AliRunLoader.h"
271 #include "AliTPCtrack.h"
272 #include "TTree.h"
273 #include "TBranch.h"
274 #include "TH2D.h"
275 #include "TCanvas.h"
276 #include "AliMagF.h"
277
278
279
280
281 void AliHBTTrackPoints::testesd(Int_t entr,const char* fname )
282 {
283   delete gAlice;
284   gAlice = 0x0;
285   AliRunLoader* rl = AliRunLoader::Open();
286   rl->LoadgAlice();
287   
288   Float_t mf = rl->GetAliRun()->Field()->SolenoidField();
289  
290   
291   TFile* fFile = TFile::Open(fname);
292   
293   if (fFile == 0x0)
294    {
295      printf("testesd: There is no suche a ESD file\n");
296      return;
297    }
298   AliESD* esd = dynamic_cast<AliESD*>(fFile->Get("0"));
299   AliESDtrack *t = esd->GetTrack(entr);
300   if (t == 0x0)
301    {
302      ::Error("testesd","Can not get track %d",entr);
303      return;
304    }
305
306   
307   Int_t N = 170;
308   AliHBTTrackPoints* tp = new AliHBTTrackPoints(N,t,mf,1.);
309   
310   Float_t xmin = -250;
311   Float_t xmax =  250;
312   
313   Float_t ymin = -250;
314   Float_t ymax =  250;
315
316   Float_t zmin = -250;
317   Float_t zmax =  250;
318   
319   TH2D* hxy = new TH2D("hxy","hxy",1000,xmin,xmax,1000,ymin,ymax);
320   TH2D* hxyt = new TH2D("hxyt","hxyt",1000,xmin,xmax,1000,ymin,ymax);
321   TH2D* hxyTR = new TH2D("hxyTR","hxyTR",1000,xmin,xmax,1000,ymin,ymax);
322
323   TH2D* hxz = new TH2D("hxz","hxz",1000,xmin,xmax,1000,zmin,zmax);
324   TH2D* hxzt = new TH2D("hxzt","hxzt",1000,xmin,xmax,1000,zmin,zmax);
325   TH2D* hxzTR = new TH2D("hxzTR","hxzTR",1000,xmin,xmax,1000,zmin,zmax);
326
327   hxyt->SetDirectory(0x0);   
328   hxy->SetDirectory(0x0);   
329   hxyTR->SetDirectory(0x0);   
330   
331   hxzt->SetDirectory(0x0);   
332   hxz->SetDirectory(0x0);   
333   hxzTR->SetDirectory(0x0);   
334
335   Float_t x,y,z;
336    
337   for (Int_t i = 0;i<N;i++)
338    {  
339       Double_t r = 84.1+i;
340       tp->PositionAt(i,x,y,z);
341       hxy->Fill(x,y);
342       hxz->Fill(x,z);
343       printf("Rdemanded %f\n",r);
344       printf("tpx %f tpy %f tpz %f Rt =%f\n", x,y,z,TMath::Hypot(x,y));
345       
346    }
347
348   rl->LoadTrackRefs();
349   TTree* treeTR = rl->TreeTR();
350   TBranch* b = treeTR->GetBranch("TPC");
351   
352   TClonesArray* trackrefs = new TClonesArray("AliTrackReference", 100);
353   AliTrackReference* tref;
354   b->SetAddress(&trackrefs);
355   
356   Int_t tlab = TMath::Abs(t->GetLabel());
357
358   Int_t netr = (Int_t)treeTR->GetEntries();
359   printf("Found %d entries in TR tree\n",netr);
360   
361   for (Int_t e = 0; e < netr; e++)
362    { 
363      treeTR->GetEntry(e);
364      tref = (AliTrackReference*)trackrefs->At(0);
365      if (tref == 0x0) continue;
366      if (tref->GetTrack() != tlab) continue;
367     
368      printf("Found %d entries in TR array\n",trackrefs->GetEntries());
369      
370      for (Int_t i = 0; i < trackrefs->GetEntries(); i++)
371       {
372         tref = (AliTrackReference*)trackrefs->At(i);
373         if (tref->GetTrack() != tlab) continue;
374         x = tref->X();
375         y = tref->Y();
376         z = tref->Z();
377         printf("Track Ref: x %f y %f z %f\n",tref->X(),tref->Y(),tref->Z());
378
379         hxzTR->Fill(x,z);
380         hxyTR->Fill(x,y);
381         for (Int_t j = 1; j < 10; j++)
382          {
383             hxyTR->Fill(x, y+j*0.1);
384             hxyTR->Fill(x, y-j*0.1);
385             hxyTR->Fill(x+j*0.1,y);
386             hxyTR->Fill(x-j*0.1,y);
387
388             hxzTR->Fill(x,z-j*0.1);
389             hxzTR->Fill(x,z+j*0.1);
390             hxzTR->Fill(x-j*0.1,z);
391             hxzTR->Fill(x+j*0.1,z);
392          }
393       }
394       break; 
395     }  
396   hxy->Draw("");
397 //  hxzt->Draw("same");
398   hxyTR->Draw("same");
399   
400   delete rl;
401 }
402
403 /***************************************************************/
404 /***************************************************************/
405 /***************************************************************/
406
407 void AliHBTTrackPoints::testtpc(Int_t entr)
408 {
409   delete gAlice;
410   gAlice = 0x0;
411   AliRunLoader* rl = AliRunLoader::Open();
412   AliLoader* l = rl->GetLoader("TPCLoader");
413   rl->LoadgAlice();
414   AliKalmanTrack::SetConvConst(100/0.299792458/0.2/rl->GetAliRun()->Field()->Factor());
415   l->LoadTracks();
416   AliTPCtrack* t = new AliTPCtrack();
417   TBranch* b=l->TreeT()->GetBranch("tracks");
418   b->SetAddress(&t);
419   l->TreeT()->GetEntry(entr);  
420   Int_t N = 160;
421   AliHBTTrackPoints* tp = new AliHBTTrackPoints(N,t,1.);
422   
423   Float_t xmin = -250;
424   Float_t xmax =  250;
425   
426   Float_t ymin = -250;
427   Float_t ymax =  250;
428
429   Float_t zmin = -250;
430   Float_t zmax =  250;
431   
432   TH2D* hxy = new TH2D("hxy","hxy",1000,xmin,xmax,1000,ymin,ymax);
433   TH2D* hxyt = new TH2D("hxyt","hxyt",1000,xmin,xmax,1000,ymin,ymax);
434   TH2D* hxyTR = new TH2D("hxyTR","hxyTR",1000,xmin,xmax,1000,ymin,ymax);
435
436   TH2D* hxz = new TH2D("hxz","hxz",1000,xmin,xmax,1000,zmin,zmax);
437   TH2D* hxzt = new TH2D("hxzt","hxzt",1000,xmin,xmax,1000,zmin,zmax);
438   TH2D* hxzTR = new TH2D("hxzTR","hxzTR",1000,xmin,xmax,1000,zmin,zmax);
439
440   hxyt->SetDirectory(0x0);   
441   hxy->SetDirectory(0x0);   
442   hxyTR->SetDirectory(0x0);   
443   
444   hxzt->SetDirectory(0x0);   
445   hxz->SetDirectory(0x0);   
446   hxzTR->SetDirectory(0x0);   
447
448   Float_t x,y,z;
449    
450   for (Int_t i = 0;i<N;i++)
451    {  
452       Double_t r = 84.1+i;
453       tp->PositionAt(i,x,y,z);
454       hxy->Fill(x,y);
455       hxz->Fill(x,z);
456       printf("Rdemanded %f\n",r);
457       printf("tpx %f tpy %f tpz %f Rt =%f\n", x,y,z,TMath::Hypot(x,y));
458       
459       //BUT they are local!!!!
460       t->PropagateTo(r);
461 //      Double_t phi = t->Phi();
462       Double_t rl = TMath::Hypot(t->GetX(),t->GetY());//real radius 
463       
464       Double_t alpha = t->GetAlpha();
465       Double_t salpha = TMath::Sin(alpha);
466       Double_t calpha = TMath::Cos(alpha);
467       x = t->GetX()*calpha - t->GetY()*salpha;
468       y = t->GetX()*salpha + t->GetY()*calpha;
469       z = t->GetZ();
470       
471       printf("tx %f ty %f tz %f Rt = %f R from XY %f\n",x,y,z,TMath::Hypot(x,y),rl);
472       
473       printf("tpz - tz %f\n",z-t->GetZ());
474       printf("\n");
475       hxyt->Fill(x,y);
476       hxzt->Fill(x,z);
477       
478    }
479
480   rl->LoadTrackRefs();
481   TTree* treeTR = rl->TreeTR();
482   b = treeTR->GetBranch("TPC");
483   
484   TClonesArray* trackrefs = new TClonesArray("AliTrackReference", 100);
485   AliTrackReference* tref;
486   b->SetAddress(&trackrefs);
487   
488   Int_t tlab = TMath::Abs(t->GetLabel());
489
490   Int_t netr = (Int_t)treeTR->GetEntries();
491   printf("Found %d entries in TR tree\n",netr);
492   
493   for (Int_t e = 0; e < netr; e++)
494    { 
495      treeTR->GetEntry(e);
496      tref = (AliTrackReference*)trackrefs->At(0);
497      if (tref == 0x0) continue;
498      if (tref->GetTrack() != tlab) continue;
499     
500      printf("Found %d entries in TR array\n",trackrefs->GetEntries());
501      
502      for (Int_t i = 0; i < trackrefs->GetEntries(); i++)
503       {
504         tref = (AliTrackReference*)trackrefs->At(i);
505         if (tref->GetTrack() != tlab) continue;
506         x = tref->X();
507         y = tref->Y();
508         z = tref->Z();
509         printf("Track Ref: x %f y %f z %f\n",tref->X(),tref->Y(),tref->Z());
510
511         hxzTR->Fill(x,z);
512         hxyTR->Fill(x,y);
513         for (Int_t j = 1; j < 10; j++)
514          {
515             hxyTR->Fill(x, y+j*0.1);
516             hxyTR->Fill(x, y-j*0.1);
517             hxyTR->Fill(x+j*0.1,y);
518             hxyTR->Fill(x-j*0.1,y);
519
520             hxzTR->Fill(x,z-j*0.1);
521             hxzTR->Fill(x,z+j*0.1);
522             hxzTR->Fill(x-j*0.1,z);
523             hxzTR->Fill(x+j*0.1,z);
524          }
525       }
526       break; 
527     }  
528   hxz->Draw("");
529 //  hxzt->Draw("same");
530   hxzTR->Draw("same");
531   
532   delete rl;
533 }