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