1 /**************************************************************************
2 * Copyright(c) 1998-2003, 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 **************************************************************************/
16 ////////////////////////////////////////////////////
17 // Stand alone tracker class //
18 // Origin: Elisabetta Crescio //
19 // e-mail: crescio@to.infn.it //
20 // tracks are saved as AliITStrackV2 objects //
21 ////////////////////////////////////////////////////
27 #include <TObjArray.h>
30 #include "AliITSclusterTable.h"
31 #include "AliITSclusterV2.h"
32 #include "AliITSgeom.h"
33 #include "AliITSRiemannFit.h"
34 #include "AliITStrackerSA.h"
35 #include "AliITStrackSA.h"
36 #include "AliITSVertexer.h"
37 #include "AliESDVertex.h"
39 #include "AliESDtrack.h"
41 ClassImp(AliITStrackerSA)
43 //____________________________________________________________________________
44 AliITStrackerSA::AliITStrackerSA():AliITStrackerV2(){
45 // Default constructor
48 //____________________________________________________________________________
49 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom):AliITStrackerV2(geom)
51 // Standard constructor (Vertex is known and passed to this obj.)
57 //____________________________________________________________________________
58 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliESDVertex *vert):AliITStrackerV2(geom)
60 // Standard constructor (Vertex is known and passed to this obj.)
66 //______________________________________________________________________
67 AliITStrackerSA::AliITStrackerSA(const AliITStrackerSA &trkr) :
68 AliITStrackerV2(trkr) {
70 // Copies are not allowed. The method is protected to avoid misuse.
71 Error("AliITStrackerSA","Copy constructor not allowed\n");
74 //______________________________________________________________________
75 AliITStrackerSA& AliITStrackerSA::operator=(const
76 AliITStrackerSA& /* trkr */){
77 // Assignment operator
78 // Assignment is not allowed. The method is protected to avoid misuse.
79 Error("= operator","Assignment operator not allowed\n");
83 //____________________________________________________________________________
84 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertexer *vertexer):AliITStrackerV2(geom)
86 // Standard constructor (Vertex is unknown - vertexer is passed to this obj)
92 //____________________________________________________________________________
93 AliITStrackerSA::AliITStrackerSA(AliITStrackerSA& tracker):AliITStrackerV2(){
95 fPhiEstimate = tracker.fPhiEstimate;
96 for(Int_t i=0;i<2;i++){
97 fPoint1[i]=tracker.fPoint1[i];
98 fPoint2[i]=tracker.fPoint2[i];
99 fPoint3[i]=tracker.fPoint3[i];
100 fPointc[i]=tracker.fPointc[i];
102 fLambdac = tracker.fLambdac;
103 fPhic = tracker.fPhic;
104 fCoef1 = tracker.fCoef1;
105 fCoef2 = tracker.fCoef2;
106 fCoef3 = tracker.fCoef3;
107 fNloop = tracker.fNloop;
108 fPhiWin = tracker.fPhiWin;
109 fLambdaWin = tracker.fLambdaWin;
110 if(tracker.fVertexer && tracker.fVert){
111 fVert = new AliESDVertex(*tracker.fVert);
114 fVert = tracker.fVert;
116 fVertexer = tracker.fVertexer;
117 fGeom = tracker.fGeom;
118 fFlagLoad = tracker.fFlagLoad;
119 fTable = tracker.fTable;
122 //____________________________________________________________________________
123 AliITStrackerSA::~AliITStrackerSA(){
125 // if fVertexer is not null, the AliESDVertex obj. is owned by this class
126 // and is deleted here
128 if(fVert)delete fVert;
134 if(fPhiWin)delete []fPhiWin;
135 if(fLambdaWin)delete []fLambdaWin;
139 //____________________________________________________________________________
140 void AliITStrackerSA::Init(){
141 // Reset all data members
143 for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
160 //_______________________________________________________________________
161 void AliITStrackerSA::ResetForFinding(){
162 // Reset data members used in all loops during track finding
164 for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
173 //____________________________________________________________________________
174 void AliITStrackerSA::FindTracks(TTree *out,Int_t evnumber){
176 /**************************************************************************
177 * This function finds primary tracks.
180 * Example: to execute function with only the ITS (no combined tracking *
181 * with TPC+ITS) and requiring 5/6 points to define a good track *
182 * call SetSixPoinbts(kFALSE) in advance and then *
183 * use: FindTracks(treein,treeout,evnumber) *
184 * to execute combined tracking, before using FindTracks, use *
186 *************************************************************************/
189 Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
194 if(fVert)delete fVert;
195 fVert = fVertexer->FindVertexForCurrentEvent(evnumber);
198 gAlice->GetEvent(evnumber);
200 Fatal("FindTracks","Vertex is missing\n");
204 Double_t primaryVertex[3];
205 Double_t errorsprimvert[3];
206 fVert->GetXYZ(primaryVertex);
207 fVert->GetSigmaXYZ(errorsprimvert);
208 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
209 Warning("FindTracks","Set errors on vertex positions x and y at 0.0001");
210 errorsprimvert[0]=0.0001;
211 errorsprimvert[1]=0.0001;
213 fVert->PrintStatus();
216 //Fill array with cluster indices for each module
218 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
219 fTable->FillArray(fITSclusters);
220 fTable->FillArrayCoorAngles();
224 //Fill tree for found tracks
225 AliITStrackV2* outrack=0;
226 TBranch* branch=out->Branch("tracks","AliITStrackV2",&outrack,32000,0);
227 if (!branch) out->Branch("tracks","AliITStrackV2",&outrack,32000,3);
228 else branch->SetAddress(&outrack);
231 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
232 for(Int_t i=0;i<fGeom->GetNlayers();i++){
233 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
235 // firstmod [i] number of the first module in the ITS layer i.
237 AliITSlayer &layer=fgLayers[0]; // first layer
239 Int_t dim=layer.GetNumberOfClusters();
240 //loop on the different windows
241 for(Int_t nloop=0;nloop<fNloop;nloop++){
242 for(Int_t ncl=0;ncl<dim;ncl++){
243 //loop starting from layer 0
246 AliITSclusterV2* cl = layer.GetCluster(ncl);
247 if(cl->IsUsed()==1) continue;
248 fPhic = fTable->GetPhiCluster(0,ncl);
249 fLambdac = fTable->GetLambdaCluster(0,ncl);
250 fPhiEstimate = fPhic;
251 AliITStrackSA* trs = new AliITStrackSA();
252 fPoint1[0]=primaryVertex[0];
253 fPoint1[1]=primaryVertex[1];
254 fPoint2[0]=fTable->GetXCluster(0,ncl);
255 fPoint2[1]=fTable->GetYCluster(0,ncl);
257 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
258 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
259 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
260 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
264 fPoint3[0] = fPointc[0];
265 fPoint3[1] = fPointc[1];
267 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
268 if(nn[1]==0 && nn[2]==0) pflag=0;
269 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
270 if(nn[2]!=0 && nn[1]==0){
272 fPoint3[0]=fPointc[0];
273 fPoint3[1]=fPointc[1];
276 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
278 if(nn[3]!=0) UpdatePoints();
279 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
281 if(nn[4]!=0) UpdatePoints();
282 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
286 Int_t numberofpoints;
287 if(fSixPoints) numberofpoints=6; //check of the candidate track
288 else numberofpoints=5; //if track is good (with the required number
289 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
290 if(nn[nnp]!=0) layOK+=1;
292 if(layOK>=numberofpoints){
293 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
295 Int_t nct = trs->GetNumberOfClustersSA();
297 Int_t index = trs->GetClusterIndexSA(nct);
298 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
299 if(kl->IsUsed()==1) kl->Use();
306 Int_t nct = tr2->GetNumberOfClusters();
309 Int_t index = tr2->GetClusterIndex(nct);
310 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
311 if(kl->IsUsed()==0) kl->Use();
316 Int_t nct = trs->GetNumberOfClustersSA();
318 Int_t index = trs->GetClusterIndexSA(nct);
319 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
320 if(kl->IsUsed()==1) kl->Use();
325 }//end loop on clusters of layer1
329 //if 5/6 points are required, second loop starting
330 //from second layer, to find tracks with point of
334 // counter for clusters on each layer
335 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
336 for(Int_t nloop=0;nloop<fNloop;nloop++){
337 AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
338 Int_t ncl2=layer2.GetNumberOfClusters();
339 while(ncl2--){ //loop starting from layer 2
342 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
343 if(cl->IsUsed()==1) continue;
344 fPhic = fTable->GetPhiCluster(1,ncl2);
345 fLambdac = fTable->GetLambdaCluster(1,ncl2);
346 fPhiEstimate = fPhic;
347 AliITStrackSA* trs = new AliITStrackSA();
349 fPoint1[0]=primaryVertex[0];
350 fPoint1[1]=primaryVertex[1];
351 fPoint2[0]=fTable->GetXCluster(1,ncl2);;
352 fPoint2[1]=fTable->GetYCluster(1,ncl2);;
354 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
355 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
356 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
357 trs,primaryVertex[2],pflag,fTable);
365 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
366 if(nn[nnp]!=0) fl+=1;
369 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
371 Int_t nct = trs->GetNumberOfClustersSA();
373 Int_t index = trs->GetClusterIndexSA(nct);
374 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
375 if(kl->IsUsed()==1) kl->Use();
381 Int_t nct = tr2->GetNumberOfClusters();
383 Int_t index = tr2->GetClusterIndex(nct);
384 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
386 if(kl->IsUsed()==0) kl->Use();
390 Int_t nct = trs->GetNumberOfClustersSA();
392 Int_t index = trs->GetClusterIndexSA(nct);
393 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
395 if(kl->IsUsed()==1) kl->Use();
399 }//end loop on clusters of layer2
402 } // if(!fSixPoints....
409 //______________________________________________________________________
410 Int_t AliITStrackerSA::FindTracks(AliESD* event){
412 // Track finder using the ESD object
415 Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
420 Double_t errorsprimvert[3];
421 Double_t primaryVertex[3];
422 event->GetVertex()->GetXYZ(primaryVertex);
423 event->GetVertex()->GetSigmaXYZ(errorsprimvert);
425 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
426 // Warning("FindTracks","Set errors on vertex positions x and y at 0.005");
427 errorsprimvert[0]=0.005;
428 errorsprimvert[1]=0.005;
431 //Fill array with cluster indices for each module
432 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
433 fTable->FillArray(fITSclusters);
434 fTable->FillArrayCoorAngles();
437 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
438 for(Int_t i=0;i<fGeom->GetNlayers();i++){
439 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
441 // firstmod [i] number of the first module in the ITS layer i.
444 AliITSlayer &layer=fgLayers[0]; // first layer
447 Int_t dim=layer.GetNumberOfClusters();
448 //loop on the different windows
449 for(Int_t nloop=0;nloop<fNloop;nloop++){
450 for(Int_t ncl=0;ncl<dim;ncl++){ //loop starting from layer 0
453 AliITSclusterV2* cl = (AliITSclusterV2*)layer.GetCluster(ncl);
455 if(cl->IsUsed()==1) continue;
456 fPhic = fTable->GetPhiCluster(0,ncl);
457 fLambdac = fTable->GetLambdaCluster(0,ncl);
458 fPhiEstimate = fPhic;
459 AliITStrackSA* trs = new AliITStrackSA();
460 fPoint1[0]=primaryVertex[0];
461 fPoint1[1]=primaryVertex[1];
462 fPoint2[0]=fTable->GetXCluster(0,ncl);
463 fPoint2[1]=fTable->GetYCluster(0,ncl);
464 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
465 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
466 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
468 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
471 fPoint3[0] = fPointc[0];
472 fPoint3[1] = fPointc[1];
474 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
475 if(nn[1]==0 && nn[2]==0) pflag=0;
476 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
477 if(nn[2]!=0 && nn[1]==0){
479 fPoint3[0]=fPointc[0];
480 fPoint3[1]=fPointc[1];
483 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
485 if(nn[3]!=0) UpdatePoints();
486 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
488 if(nn[4]!=0) UpdatePoints();
489 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
493 Int_t numberofpoints;
494 if(fSixPoints) numberofpoints=6; //check of the candidate track
495 else numberofpoints=5; //if track is good (with the required number
496 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
497 if(nn[nnp]!=0) layOK+=1;
499 if(layOK>=numberofpoints){
500 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
502 Int_t nct = trs->GetNumberOfClustersSA();
504 Int_t index = trs->GetClusterIndexSA(nct);
505 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
506 if(kl->IsUsed()==1) kl->Use();
511 AliESDtrack outtrack;
512 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
513 event->AddTrack(&outtrack);
516 Int_t nct = tr2->GetNumberOfClusters();
518 Int_t index = tr2->GetClusterIndex(nct);
519 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
520 if(kl->IsUsed()==0) kl->Use();
524 Int_t nct = trs->GetNumberOfClustersSA();
526 Int_t index = trs->GetClusterIndexSA(nct);
527 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
528 if(kl->IsUsed()==1) kl->Use();
534 }//end loop on clusters of layer1
538 //if 5/6 points are required, second loop starting
539 //from second layer, to find tracks with point of
543 // counter for clusters on each layer
544 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
545 for(Int_t nloop=0;nloop<fNloop;nloop++){
546 AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
547 Int_t ncl2=layer2.GetNumberOfClusters();
548 while(ncl2--){ //loop starting from layer 2
551 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
552 if(cl->IsUsed()==1) continue;
553 fPhic = fTable->GetPhiCluster(1,ncl2);
554 fLambdac = fTable->GetLambdaCluster(1,ncl2);
555 fPhiEstimate = fPhic;
556 AliITStrackSA* trs = new AliITStrackSA();
557 fPoint1[0]=primaryVertex[0];
558 fPoint1[1]=primaryVertex[1];
559 fPoint2[0]=fTable->GetXCluster(1,ncl2);
560 fPoint2[1]=fTable->GetYCluster(1,ncl2);
562 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
563 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
564 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
565 trs,primaryVertex[2],pflag,fTable);
573 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
574 if(nn[nnp]!=0) fl+=1;
577 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
580 Int_t nct = trs->GetNumberOfClustersSA();
582 Int_t index = trs->GetClusterIndexSA(nct);
583 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
584 if(kl->IsUsed()==1) kl->Use();
589 AliESDtrack outtrack;
590 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
591 event->AddTrack(&outtrack);
593 Int_t nct = tr2->GetNumberOfClusters();
595 Int_t index = tr2->GetClusterIndex(nct);
596 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
598 if(kl->IsUsed()==0) kl->Use();
602 Int_t nct = trs->GetNumberOfClustersSA();
604 Int_t index = trs->GetClusterIndexSA(nct);
605 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
607 if(kl->IsUsed()==1) kl->Use();
611 }//end loop on clusters of layer2
618 Info("FindTracks","Number of found tracks: %d",event->GetNumberOfTracks());
625 //________________________________________________________________________
627 AliITStrackV2* AliITStrackerSA::FitTrack(AliITStrackSA* tr,Double_t *primaryVertex,Double_t *errorsprimvert){
628 //fit of the found track
631 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
632 for(Int_t i=0;i<fGeom->GetNlayers();i++){
633 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
635 TObjArray* listoftracks=new TObjArray(0,0);
636 AliITStrackV2* otrack2;
637 Int_t nclusters = tr->GetNumberOfClustersSA();
638 TObjArray** listlayer = new TObjArray*[fGeom->GetNlayers()];
639 for(Int_t i=0;i<fGeom->GetNlayers();i++){
640 listlayer[i] = new TObjArray(0,0);
650 Int_t * nnn = new Int_t[fGeom->GetNlayers()];
651 for(Int_t i=0;i<fGeom->GetNlayers();i++)nnn[i]=0;
653 for(Int_t ncl=0;ncl<nclusters;ncl++){
654 Int_t index = tr->GetClusterIndexSA(ncl);
655 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(index);
657 if(cl->IsUsed()==1) cl->Use();
658 Int_t lay = (index & 0xf0000000) >> 28;
659 if(lay==0) { listlayer[0]->AddLast(cl); clind0[nnn[0]]=index;nnn[0]++;}
660 if(lay==1) { listlayer[1]->AddLast(cl); clind1[nnn[1]]=index;nnn[1]++;}
661 if(lay==2) { listlayer[2]->AddLast(cl); clind2[nnn[2]]=index;nnn[2]++;}
662 if(lay==3) { listlayer[3]->AddLast(cl); clind3[nnn[3]]=index;nnn[3]++;}
663 if(lay==4) { listlayer[4]->AddLast(cl); clind4[nnn[4]]=index;nnn[4]++;}
664 if(lay==5) { listlayer[5]->AddLast(cl); clind5[nnn[5]]=index;nnn[5]++;}
669 Int_t * end = new Int_t[fGeom->GetNlayers()];
670 for(Int_t i=0;i<fGeom->GetNlayers();i++){
671 if(listlayer[i]->GetEntries()==0) end[i]=1;
672 else end[i]=listlayer[i]->GetEntries();
675 for(Int_t l1=0;l1<end[0];l1++){//loop on layer 1
676 AliITSclusterV2* cl0 = (AliITSclusterV2*)listlayer[0]->At(l1);
677 TVector3** recp = new TVector3*[3];
678 TVector3** errs = new TVector3*[3];
679 recp[0] = new TVector3(primaryVertex[0],primaryVertex[1],primaryVertex[2]);
680 errs[0] = new TVector3(errorsprimvert[0],errorsprimvert[1],errorsprimvert[2]);
681 Double_t x1,y1,z1,sx1,sy1,sz1;
682 Double_t x2,y2,z2,sx2,sy2,sz2;
683 AliITSclusterV2* p1=0;
684 AliITSclusterV2* p2=0;
685 Int_t index1=clind0[l1];
687 for(Int_t l2=0;l2<end[1];l2++){//loop on layer 2
688 AliITSclusterV2* cl1 = (AliITSclusterV2*)listlayer[1]->At(l2);
690 for(Int_t l3=0;l3<end[2];l3++){ //loop on layer 3
691 AliITSclusterV2* cl2 = (AliITSclusterV2*)listlayer[2]->At(l3);
693 if(cl0==0 && cl1!=0) {
694 p2 = cl2;index1=clind2[l3];
698 if(cl0!=0 && cl1==0){
700 p2=cl2;index2=clind2[l3];
702 if(cl0!=0 && cl1!=0){
706 Int_t lay1=(index1 & 0xf0000000) >> 28;
707 Int_t cln1=(index1 & 0x0fffffff) >> 00;
708 Int_t lay2=(index2 & 0xf0000000) >> 28;
709 Int_t cln2=(index2 & 0x0fffffff) >> 00;
710 x1 = fTable->GetXCluster(lay1,cln1);
711 x2 = fTable->GetXCluster(lay2,cln2);
712 y1 = fTable->GetYCluster(lay1,cln1);
713 y2 = fTable->GetYCluster(lay2,cln2);
714 z1 = fTable->GetZCluster(lay1,cln1);
715 z2 = fTable->GetZCluster(lay2,cln2);
716 sx1 = fTable->GetXClusterError(lay1,cln1);
717 sx2 = fTable->GetXClusterError(lay2,cln2);
718 sy1 = fTable->GetYClusterError(lay1,cln1);
719 sy2 = fTable->GetYClusterError(lay2,cln2);
720 sz1 = fTable->GetZClusterError(lay1,cln1);
721 sz2 = fTable->GetZClusterError(lay2,cln2);
722 Double_t phi1 = fTable->GetPhiCluster(lay1,cln1);
723 Int_t module1 = p1->GetDetectorIndex()+firstmod[0];
724 recp[1] = new TVector3(x1,y1,z1);
725 errs[1] = new TVector3(sx1,sy1,sz1);
726 recp[2] = new TVector3(x2,y2,z2);
727 errs[2] = new TVector3(sx2,sy2,sz2);
729 //fit on the Riemann sphere
730 Float_t seed1,seed2,seed3;
731 AliITSRiemannFit fit;
732 Int_t rf = fit.FitHelix(3,recp,errs,seed1,seed2,seed3); //this gives phi,tgl,curvature to start Kalman Filter
738 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
739 phi=seed1+1.5*TMath::Pi();
742 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
743 phi=seed1+0.5*TMath::Pi();
746 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
747 phi=seed1-0.5*TMath::Pi();
752 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
753 phi=seed1+0.5*TMath::Pi();
756 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
757 phi=seed1-0.5*TMath::Pi();
760 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
761 phi=seed1-1.5*TMath::Pi();
766 Int_t layer,ladder,detector;
767 fGeom->GetModuleId(module1,layer,ladder,detector);
768 Float_t yclu1 = p1->GetY();
769 Float_t zclu1 = p1->GetZ();
770 Double_t cv=Curvature(primaryVertex[0],primaryVertex[1],x1,y1,x2,y2);
772 for(Int_t l4=0;l4<end[3];l4++){ //loop on layer 4
773 AliITSclusterV2* cl3 = (AliITSclusterV2*)listlayer[3]->At(l4);
774 for(Int_t l5=0;l5<end[4];l5++){ //loop on layer 5
775 AliITSclusterV2* cl4 = (AliITSclusterV2*)listlayer[4]->At(l5);
776 for(Int_t l6=0;l6<end[5];l6++){ //loop on layer 6
777 AliITSclusterV2* cl5 = (AliITSclusterV2*)listlayer[5]->At(l6);
778 AliITStrackSA* trac = new AliITStrackSA(layer,ladder,detector,yclu1,zclu1,phi,tgl,cv,1);
780 if(cl5!=0) trac->AddClusterV2(5,(clind5[l6] & 0x0fffffff)>>0);
781 if(cl4!=0) trac->AddClusterV2(4,(clind4[l5] & 0x0fffffff)>>0);
782 if(cl3!=0) trac->AddClusterV2(3,(clind3[l4] & 0x0fffffff)>>0);
783 if(cl2!=0) trac->AddClusterV2(2,(clind2[l3] & 0x0fffffff)>>0);
784 if(cl1!=0) trac->AddClusterV2(1,(clind1[l2] & 0x0fffffff)>>0);
785 if(cl0!=0) trac->AddClusterV2(0,(clind0[l1] & 0x0fffffff)>>0);
787 //fit with Kalman filter using AliITStrackerV2::RefitAt()
789 AliITStrackV2* ot = new AliITStrackV2(*trac);
791 ot->ResetCovariance();
794 if(RefitAt(49.,ot,trac)){ //fit from layer 1 to layer 6
796 otrack2 = new AliITStrackV2(*ot);
797 otrack2->ResetCovariance();
798 otrack2->ResetClusters();
799 //fit from layer 6 to layer 1
800 if(RefitAt(3.7,otrack2,ot)) listoftracks->AddLast(otrack2);
810 for(Int_t i=1;i<3;i++){
824 Int_t dim=listoftracks->GetEntries();
827 for(Int_t i=0;i<fGeom->GetNlayers();i++){
834 AliITStrackV2* otrack =(AliITStrackV2*)FindTrackLowChiSquare(listoftracks,dim);
836 if(otrack==0) return 0;
837 Int_t * indexc = new Int_t[fGeom->GetNlayers()];
838 for(Int_t i=0;i<fGeom->GetNlayers();i++) indexc[i]=0;
839 for(Int_t nind=0;nind<otrack->GetNumberOfClusters();nind++){
840 indexc[nind] = otrack->GetClusterIndex(nind);
842 AliITSclusterV2* cl0 = (AliITSclusterV2*)GetCluster(indexc[0]);
843 AliITSclusterV2* cl1 = (AliITSclusterV2*)GetCluster(indexc[1]);
844 AliITSclusterV2* cl2 = (AliITSclusterV2*)GetCluster(indexc[2]);
845 AliITSclusterV2* cl3 = (AliITSclusterV2*)GetCluster(indexc[3]);
846 AliITSclusterV2* cl4 = (AliITSclusterV2*)GetCluster(indexc[4]);
847 Int_t labl[3]={-1,-1,-1};
848 if(otrack->GetNumberOfClusters()==fGeom->GetNlayers()){
849 AliITSclusterV2* cl5 = (AliITSclusterV2*)GetCluster(indexc[5]);
850 labl[0]=cl5->GetLabel(0);
851 labl[1]=cl5->GetLabel(1);
852 labl[2]=cl5->GetLabel(2);
855 if(otrack->GetNumberOfClusters()==(fGeom->GetNlayers()-1)){
860 Int_t numberofpoints;
861 if(fSixPoints) numberofpoints=6;
862 else numberofpoints=5;
863 Int_t label = Label(cl0->GetLabel(0),cl1->GetLabel(0),
864 cl2->GetLabel(0),cl3->GetLabel(0),
865 cl4->GetLabel(0),labl[0],
866 cl0->GetLabel(1),cl1->GetLabel(1),
867 cl2->GetLabel(1),cl3->GetLabel(1),
868 cl4->GetLabel(1),labl[1],
869 cl0->GetLabel(2),cl1->GetLabel(2),
870 cl2->GetLabel(2),cl3->GetLabel(2),
871 cl4->GetLabel(2),labl[2],numberofpoints);
873 otrack->SetLabel(label);
875 for(Int_t i=0;i<fGeom->GetNlayers();i++){
884 //_______________________________________________________________________
885 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,TTree* treev2, TTree* clustertree){
886 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
887 //(or AliITStrackV2 tracks found with function FindTracks of this class)
891 if(fVert)delete fVert;
892 fVert = fVertexer->FindVertexForCurrentEvent(evnum);
895 gAlice->GetEvent(evnum);
897 Fatal("FindTracks","Vertex is missing\n");
901 Double_t primaryVertex[3];
902 fVert->GetXYZ(primaryVertex);
905 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
906 fTable->FillArray(fITSclusters);
907 fTable->FillArrayCoorAngles();
909 SetEventNumber(evnum);
910 if(GetFlagLoadedClusters()==0) LoadClusters(clustertree);
911 SetFlagLoadedClusters(1);
913 TBranch* bra = (TBranch*)treev2->GetBranch("tracks");
914 if(!bra) Warning("UseFoundTracksV2","No branch for track tree");
915 AliITStrackV2* ttrrt = new AliITStrackV2;
916 bra->SetAddress(&ttrrt);
918 for(Int_t nj=0;nj<treev2->GetEntries();nj++){
919 treev2->GetEvent(nj);
920 Int_t ncl = ttrrt->GetNumberOfClusters();
921 for(Int_t k=0;k<ncl;k++){
922 Int_t index = ttrrt->GetClusterIndex(k);
923 AliITSclusterV2* clui = (AliITSclusterV2*)GetCluster(index);
924 if(clui->IsUsed()==0) clui->Use();
931 //_______________________________________________________________________
932 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,AliESD *event, TTree*
934 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
938 if(fVert)delete fVert;
939 fVert = fVertexer->FindVertexForCurrentEvent(evnum);
942 gAlice->GetEvent(evnum);
944 Fatal("FindTracks","Vertex is missing\n");
948 Double_t primaryVertex[3];
949 fVert->GetXYZ(primaryVertex);
952 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
953 fTable->FillArray(fITSclusters);
954 fTable->FillArrayCoorAngles();
956 SetEventNumber(evnum);
957 if(GetFlagLoadedClusters()==0) LoadClusters(clustertree);
958 SetFlagLoadedClusters(1);
959 Int_t ntracks = event->GetNumberOfTracks();
961 AliESDtrack *esd=event->GetTrack(ntracks);
962 if ((esd->GetStatus()&
963 AliESDtrack::kITSin|AliESDtrack::kTPCin)==0) continue;
965 Int_t ncl = esd->GetITSclusters(idx);
966 for(Int_t clu=0; clu<ncl; clu++){
967 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(idx[clu]);
968 if(cl->IsUsed()==0) cl->Use();
972 Info("UseFoundTracksV2","Clusters of tracks prolonged from TPC deleted");
977 //_______________________________________________________
978 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Int_t pflag,AliITSclusterTable* table){
979 //function used to to find the clusters associated to the track
981 AliITSlayer &lay = fgLayers[layer];
982 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
983 for(Int_t i=0;i<fGeom->GetNlayers();i++){
984 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
988 Float_t cx1,cx2,cy1,cy2;
989 FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
990 Int_t fun = FindIntersection(fCoef1,fCoef2,fCoef3,-(lay.GetR()*lay.GetR()),cx1,cy1,cx2,cy2);
995 Double_t fi1 =TMath::ATan2(cy1,cx1);
996 Double_t fi2 =TMath::ATan2(cy2,cx2);
997 fPhiEstimate = ChoosePoint(fi1,fi2,fPhic);
1000 Double_t zed = TMath::Tan(fLambdac)*lay.GetR()+zvertex;
1001 Double_t zed1 = TMath::Tan(fLambdac+lambdawindow)*lay.GetR()+zvertex;
1002 Double_t zed2 = TMath::Tan(fLambdac-lambdawindow)*lay.GetR()+zvertex;
1004 Double_t fi = fPhiEstimate;
1005 Int_t nmod = lay.FindDetectorIndex(fi,zed);
1010 nmod += firstmod[layer];
1012 Int_t nm[8]={0,0,0,0,0,0,0,0};
1013 nm[0] = lay.FindDetectorIndex(fi+phiwindow,zed);
1014 nm[1] = lay.FindDetectorIndex(fi-phiwindow,zed);
1015 nm[2] = lay.FindDetectorIndex(fi,zed1);
1016 nm[3] = lay.FindDetectorIndex(fi,zed2);
1017 nm[4] = lay.FindDetectorIndex(fi+phiwindow,zed1);
1018 nm[5] = lay.FindDetectorIndex(fi-phiwindow,zed1);
1019 nm[6] = lay.FindDetectorIndex(fi+phiwindow,zed2);
1020 nm[7] = lay.FindDetectorIndex(fi-phiwindow,zed2);
1023 TArrayI* array =(TArrayI*)table->GetListOfClusters(nmod);
1024 TArrayI* list = new TArrayI(array->GetSize());
1025 for(Int_t i=0;i<array->GetSize();i++){
1026 Int_t in=(Int_t)array->At(i);
1031 // Info("SearchClusters", "layer %d, module %d", layer, nmod);
1032 for(Int_t ii=0;ii<8;ii++){
1033 if(nm[ii]!=nmod && nm[ii]>=0){
1034 TArrayI* ar =(TArrayI*)table->GetListOfClusters(nm[ii]+firstmod[layer]);
1035 list->Set(list->GetSize()+ar->GetSize());
1036 for(Int_t j=0;j<ar->GetSize();j++){
1037 Int_t in=(Int_t)ar->At(j);
1043 for(Int_t i=0;i<list->GetSize();i++){
1044 Int_t index = (Int_t)list->At(i);
1045 AliITSclusterV2* cllay = lay.GetCluster(index);
1046 if(cllay==0) continue;
1047 if(cllay->IsUsed()==1) continue;
1048 Double_t phi = fTable->GetPhiCluster(layer,index);
1049 Double_t lambda= fTable->GetLambdaCluster(layer,index);
1051 if(TMath::Abs(fLambdac-lambda)<lambdawindow &&
1052 TMath::Abs(fPhiEstimate-phi)<phiwindow){
1055 if(trs->GetNumberOfClustersSA()==20){
1060 trs->AddClusterSA(layer,index);
1063 fPointc[0]=fTable->GetXCluster(layer,index);
1064 fPointc[1]=fTable->GetYCluster(layer,index);
1076 //________________________________________________________________
1077 void AliITStrackerSA::UpdatePoints(){
1078 //update of points for the estimation of the curvature
1080 fPoint1[0]=fPoint2[0];
1081 fPoint1[1]=fPoint2[1];
1082 fPoint2[0]=fPoint3[0];
1083 fPoint2[1]=fPoint3[1];
1084 fPoint3[0]=fPointc[0];
1085 fPoint3[1]=fPointc[1];
1092 //___________________________________________________________________
1093 Int_t AliITStrackerSA::FindEquation(Float_t x1, Float_t y1, Float_t x2, Float_t y2, Float_t x3, Float_t y3,Float_t& a, Float_t& b, Float_t& c){
1095 //given (x,y) of three recpoints (in global coordinates)
1096 //returns the parameters a,b,c of circonference x*x + y*y +a*x + b*y +c
1098 Float_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1099 if(den==0) return 0;
1100 a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1101 b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1102 c = -x1*x1-y1*y1-a*x1-b*y1;
1105 //__________________________________________________________________________
1106 Int_t AliITStrackerSA::FindIntersection(Float_t a1, Float_t b1, Float_t c1, Float_t c2,Float_t& x1,Float_t& y1, Float_t& x2, Float_t& y2){
1108 //Finds the intersection between the circonference of the track and the circonference centered in (0,0) represented by one layer
1109 //c2 is -rlayer*rlayer
1113 Float_t aA = (b1*b1)/(a1*a1)+1;
1114 Float_t bB = (-2*m*b1/(a1*a1));
1115 Float_t cC = c2+(m*m)/(a1*a1);
1116 if((bB*bB-4*aA*cC)<0) return 0;
1118 y1 = (-bB+TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA);
1119 y2 = (-bB-TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA);
1120 x1 = (c2-c1-b1*y1)/a1;
1121 x2 = (c2-c1-b1*y2)/a1;
1125 //____________________________________________________________________
1126 Double_t AliITStrackerSA::Curvature(Double_t x1,Double_t y1,Double_t
1127 x2,Double_t y2,Double_t x3,Double_t y3){
1129 //calculates the curvature of track
1130 Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1131 if(den==0) return 0;
1132 Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1133 Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1134 Double_t c = -x1*x1-y1*y1-a*x1-b*y1;
1137 if((a*a+b*b-4*c)<0) return 0;
1138 Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.;
1139 if(rad==0) return 0;
1141 if((x1>0 && y1>0 && x1<xc)) rad*=-1;
1142 if((x1<0 && y1>0 && x1<xc)) rad*=-1;
1143 // if((x1<0 && y1<0 && x1<xc)) rad*=-1;
1144 // if((x1>0 && y1<0 && x1<xc)) rad*=-1;
1149 //____________________________________________________________________
1150 Double_t AliITStrackerSA::ChoosePoint(Double_t p1, Double_t p2, Double_t pp){
1152 //Returns the point closest to pp
1154 Double_t diff1 = p1-pp;
1155 Double_t diff2 = p2-pp;
1157 if(TMath::Abs(diff1)<TMath::Abs(diff2)) fPhiEstimate=p1;
1158 else fPhiEstimate=p2;
1159 return fPhiEstimate;
1164 //_________________________________________________________________
1165 AliITStrackV2* AliITStrackerSA::FindTrackLowChiSquare(TObjArray* tracklist, Int_t dim) const {
1166 // returns track with lowes chi square
1168 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(0);
1171 if(dim==0) return 0;
1172 Double_t * chi2 = new Double_t[dim];
1173 Int_t * index = new Int_t[dim];
1174 for(Int_t i=0;i<dim;i++){
1175 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(i);
1176 chi2[i]=trk->GetChi2();
1180 Int_t w=0;Double_t value;
1183 for(Int_t j=w+1;j<dim;j++){
1184 if(chi2[w]<chi2[j]){
1196 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(index[dim-1]);
1203 //__________________________________________________________
1204 Int_t AliITStrackerSA::FindLabel(Int_t l1, Int_t l2, Int_t l3, Int_t l4, Int_t l5, Int_t l6){
1206 //function used to determine the track label
1208 Int_t lb[6] = {l1,l2,l3,l4,l5,l6};
1209 Int_t aa[6]={1,1,1,1,1,1};
1212 Int_t k=0;Int_t w=0;Int_t num=6;
1213 if(lb[5]==-1) num=5;
1217 for(Int_t i=k+1;i<num;i++){
1219 if(lb[k]==lb[i] && aa[k]!=0){
1230 for(Int_t j=0;j<6;j++){
1242 if(num==6) return lb[5];
1246 //_____________________________________________________________________________
1247 Int_t AliITStrackerSA::Label(Int_t gl1, Int_t gl2, Int_t gl3, Int_t gl4, Int_t gl5, Int_t gl6,Int_t gl7, Int_t gl8, Int_t gl9, Int_t gl10,Int_t gl11,
1248 Int_t gl12, Int_t gl13, Int_t gl14,Int_t gl15, Int_t gl16, Int_t gl17, Int_t gl18, Int_t numberofpoints){
1251 //function used to assign label to the found track. If track is fake, the label is negative
1253 Int_t lb0[6] = {gl1,gl2,gl3,gl4,gl5,gl6};
1254 Int_t lb1[6] = {gl7,gl8,gl9,gl10,gl11,gl12};
1255 Int_t lb2[6] = {gl13,gl14,gl15,gl16,gl17,gl18};
1256 Int_t ll=FindLabel(lb0[0],lb0[1],lb0[2],lb0[3],lb0[4],lb0[5]);
1257 Int_t lflag=0;Int_t num=6;
1258 if(lb0[5]==-1 && lb1[5]==-1 && lb2[5]==-1) num=5;
1260 for(Int_t i=0;i<num;i++){
1261 if(lb0[i]==ll || lb1[i]==ll || lb2[i]==ll) lflag+=1;
1264 if(lflag>=numberofpoints) return ll;
1270 //_____________________________________________________________________________
1271 void AliITStrackerSA::SetWindowSizes(Int_t n, Double_t *phi, Double_t *lam){
1272 // Set sizes of the phi and lambda windows used for track finding
1274 if(phi){ // user defined values
1275 fPhiWin = new Double_t[fNloop];
1276 fLambdaWin = new Double_t[fNloop];
1277 for(Int_t k=0;k<fNloop;k++){
1279 fLambdaWin[k]=lam[k];
1282 else { // default values
1284 Double_t phid[40] = {0.001,0.0015,0.002,0.0023,0.0025,0.0027,0.003,
1285 0.0033,0.0035,0.0037,0.004,0.0043,0.0045,0.0047,
1286 0.005,0.0053,0.0055,
1287 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1288 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1289 0.01,0.0105,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1291 Double_t lambdad[40] = {0.001,0.0015,0.002,0.0023,0.0025,0.0027,0.003,
1292 0.0033,0.0035,0.0037,0.004,0.0043,0.0045,0.0047,
1293 0.005,0.0053,0.0055,
1294 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1295 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1296 0.01,0.015,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1301 Warning("SetWindowSizes","Number of loop forced to the default value %d",fNloop);
1304 fPhiWin = new Double_t[fNloop];
1305 fLambdaWin = new Double_t[fNloop];
1307 for(Int_t k=0;k<fNloop;k++){
1309 fLambdaWin[k]=lambdad[k];