2 /**************************************************************************
3 * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
17 ////////////////////////////////////////////////////
18 // Stand alone tracker class //
19 // Origin: Elisabetta Crescio //
20 // e-mail: crescio@to.infn.it //
21 // tracks are saved as AliITStrackV2 objects //
22 ////////////////////////////////////////////////////
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 fTable = tracker.fTable;
119 fListOfTracks = tracker.fListOfTracks;
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;
133 if(fPhiWin)delete []fPhiWin;
134 if(fLambdaWin)delete []fLambdaWin;
136 fListOfTracks->Delete();
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;}
158 fListOfTracks=new TObjArray(0,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;}
172 fListOfTracks->Delete();
174 //____________________________________________________________________________
175 void AliITStrackerSA::FindTracks(TTree *out,Int_t evnumber){
177 /**************************************************************************
178 * This function finds primary tracks.
181 * Example: to execute function with only the ITS (no combined tracking *
182 * with TPC+ITS) and requiring 5/6 points to define a good track *
183 * call SetSixPoinbts(kFALSE) in advance and then *
184 * use: FindTracks(treein,treeout,evnumber) *
185 * to execute combined tracking, before using FindTracks, use *
187 *************************************************************************/
190 Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
195 if(fVert)delete fVert;
196 fVert = fVertexer->FindVertexForCurrentEvent(evnumber);
199 gAlice->GetEvent(evnumber);
201 Fatal("FindTracks","Vertex is missing\n");
205 Double_t primaryVertex[3];
206 Double_t errorsprimvert[3];
207 fVert->GetXYZ(primaryVertex);
208 fVert->GetSigmaXYZ(errorsprimvert);
209 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
210 Warning("FindTracks","Set errors on vertex positions x and y at 0.0001");
211 errorsprimvert[0]=0.0001;
212 errorsprimvert[1]=0.0001;
214 fVert->PrintStatus();
217 //Fill array with cluster indices for each module
219 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
220 fTable->FillArray(fITSclusters);
221 fTable->FillArrayCoorAngles();
225 //Fill tree for found tracks
226 AliITStrackV2* outrack=0;
227 TBranch* branch=out->Branch("tracks","AliITStrackV2",&outrack,32000,0);
228 if (!branch) out->Branch("tracks","AliITStrackV2",&outrack,32000,3);
229 else branch->SetAddress(&outrack);
232 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
233 for(Int_t i=0;i<fGeom->GetNlayers();i++){
234 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
236 // firstmod [i] number of the first module in the ITS layer i.
238 AliITSlayer &layer=fgLayers[0]; // first layer
240 Int_t dim=layer.GetNumberOfClusters();
241 //loop on the different windows
242 for(Int_t nloop=0;nloop<fNloop;nloop++){
243 for(Int_t ncl=0;ncl<dim;ncl++){
244 //loop starting from layer 0
247 AliITSclusterV2* cl = layer.GetCluster(ncl);
248 if(cl->IsUsed()==1) continue;
249 fPhic = fTable->GetPhiCluster(0,ncl);
250 fLambdac = fTable->GetLambdaCluster(0,ncl);
251 fPhiEstimate = fPhic;
252 AliITStrackSA* trs = new AliITStrackSA();
253 fPoint1[0]=primaryVertex[0];
254 fPoint1[1]=primaryVertex[1];
255 fPoint2[0]=fTable->GetXCluster(0,ncl);
256 fPoint2[1]=fTable->GetYCluster(0,ncl);
258 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
259 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
260 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
261 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
265 fPoint3[0] = fPointc[0];
266 fPoint3[1] = fPointc[1];
268 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
269 if(nn[1]==0 && nn[2]==0) pflag=0;
270 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
271 if(nn[2]!=0 && nn[1]==0){
273 fPoint3[0]=fPointc[0];
274 fPoint3[1]=fPointc[1];
277 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
279 if(nn[3]!=0) UpdatePoints();
280 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
282 if(nn[4]!=0) UpdatePoints();
283 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
287 Int_t numberofpoints;
288 if(fSixPoints) numberofpoints=6; //check of the candidate track
289 else numberofpoints=5; //if track is good (with the required number
290 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
291 if(nn[nnp]!=0) layOK+=1;
293 if(layOK>=numberofpoints){
294 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
296 Int_t nct = trs->GetNumberOfClustersSA();
298 Int_t index = trs->GetClusterIndexSA(nct);
299 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
300 if(kl->IsUsed()==1) kl->Use();
307 Int_t nct = tr2->GetNumberOfClusters();
310 Int_t index = tr2->GetClusterIndex(nct);
311 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
312 if(kl->IsUsed()==0) kl->Use();
317 Int_t nct = trs->GetNumberOfClustersSA();
319 Int_t index = trs->GetClusterIndexSA(nct);
320 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
321 if(kl->IsUsed()==1) kl->Use();
326 }//end loop on clusters of layer1
330 //if 5/6 points are required, second loop starting
331 //from second layer, to find tracks with point of
335 // counter for clusters on each layer
336 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
337 for(Int_t nloop=0;nloop<fNloop;nloop++){
338 AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
339 Int_t ncl2=layer2.GetNumberOfClusters();
340 while(ncl2--){ //loop starting from layer 2
343 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
344 if(cl->IsUsed()==1) continue;
345 fPhic = fTable->GetPhiCluster(1,ncl2);
346 fLambdac = fTable->GetLambdaCluster(1,ncl2);
347 fPhiEstimate = fPhic;
348 AliITStrackSA* trs = new AliITStrackSA();
350 fPoint1[0]=primaryVertex[0];
351 fPoint1[1]=primaryVertex[1];
352 fPoint2[0]=fTable->GetXCluster(1,ncl2);;
353 fPoint2[1]=fTable->GetYCluster(1,ncl2);;
355 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
356 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
357 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
358 trs,primaryVertex[2],pflag,fTable);
366 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
367 if(nn[nnp]!=0) fl+=1;
370 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
372 Int_t nct = trs->GetNumberOfClustersSA();
374 Int_t index = trs->GetClusterIndexSA(nct);
375 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
376 if(kl->IsUsed()==1) kl->Use();
382 Int_t nct = tr2->GetNumberOfClusters();
384 Int_t index = tr2->GetClusterIndex(nct);
385 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
387 if(kl->IsUsed()==0) kl->Use();
391 Int_t nct = trs->GetNumberOfClustersSA();
393 Int_t index = trs->GetClusterIndexSA(nct);
394 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
396 if(kl->IsUsed()==1) kl->Use();
400 }//end loop on clusters of layer2
403 } // if(!fSixPoints....
410 //______________________________________________________________________
411 Int_t AliITStrackerSA::FindTracks(AliESD* event){
413 // Track finder using the ESD object
417 Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
422 Double_t errorsprimvert[3];
423 Double_t primaryVertex[3];
424 event->GetVertex()->GetXYZ(primaryVertex);
425 event->GetVertex()->GetSigmaXYZ(errorsprimvert);
427 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
428 // Warning("FindTracks","Set errors on vertex positions x and y at 0.005");
429 errorsprimvert[0]=0.005;
430 errorsprimvert[1]=0.005;
433 //Fill array with cluster indices for each module
435 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
436 fTable->FillArray(fITSclusters);
437 fTable->FillArrayCoorAngles();
440 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
441 for(Int_t i=0;i<fGeom->GetNlayers();i++){
442 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
444 // firstmod [i] number of the first module in the ITS layer i.
447 AliITSlayer &layer=fgLayers[0]; // first layer
450 Int_t dim=layer.GetNumberOfClusters();
451 //loop on the different windows
452 for(Int_t nloop=0;nloop<fNloop;nloop++){
453 for(Int_t ncl=0;ncl<dim;ncl++){ //loop starting from layer 0
456 AliITSclusterV2* cl = (AliITSclusterV2*)layer.GetCluster(ncl);
458 if(cl->IsUsed()==1) continue;
459 fPhic = fTable->GetPhiCluster(0,ncl);
460 fLambdac = fTable->GetLambdaCluster(0,ncl);
461 fPhiEstimate = fPhic;
462 AliITStrackSA* trs = new AliITStrackSA();
463 fPoint1[0]=primaryVertex[0];
464 fPoint1[1]=primaryVertex[1];
465 fPoint2[0]=fTable->GetXCluster(0,ncl);
466 fPoint2[1]=fTable->GetYCluster(0,ncl);
467 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
468 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
469 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
471 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
474 fPoint3[0] = fPointc[0];
475 fPoint3[1] = fPointc[1];
477 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
478 if(nn[1]==0 && nn[2]==0) pflag=0;
479 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
480 if(nn[2]!=0 && nn[1]==0){
482 fPoint3[0]=fPointc[0];
483 fPoint3[1]=fPointc[1];
486 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
488 if(nn[3]!=0) UpdatePoints();
489 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
491 if(nn[4]!=0) UpdatePoints();
492 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
496 Int_t numberofpoints;
497 if(fSixPoints) numberofpoints=6; //check of the candidate track
498 else numberofpoints=5; //if track is good (with the required number
499 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
500 if(nn[nnp]!=0) layOK+=1;
502 if(layOK>=numberofpoints){
503 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
505 Int_t nct = trs->GetNumberOfClustersSA();
507 Int_t index = trs->GetClusterIndexSA(nct);
508 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
509 if(kl->IsUsed()==1) kl->Use();
514 AliESDtrack outtrack;
515 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
516 event->AddTrack(&outtrack);
518 Int_t nct = tr2->GetNumberOfClusters();
520 Int_t index = tr2->GetClusterIndex(nct);
521 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
522 if(kl->IsUsed()==0) kl->Use();
527 Int_t nct = trs->GetNumberOfClustersSA();
529 Int_t index = trs->GetClusterIndexSA(nct);
530 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
531 if(kl->IsUsed()==1) kl->Use();
537 }//end loop on clusters of layer1
541 //if 5/6 points are required, second loop starting
542 //from second layer, to find tracks with point of
546 // counter for clusters on each layer
547 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
548 for(Int_t nloop=0;nloop<fNloop;nloop++){
549 AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
550 Int_t ncl2=layer2.GetNumberOfClusters();
551 while(ncl2--){ //loop starting from layer 2
554 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
555 if(cl->IsUsed()==1) continue;
556 fPhic = fTable->GetPhiCluster(1,ncl2);
557 fLambdac = fTable->GetLambdaCluster(1,ncl2);
558 fPhiEstimate = fPhic;
559 AliITStrackSA* trs = new AliITStrackSA();
560 fPoint1[0]=primaryVertex[0];
561 fPoint1[1]=primaryVertex[1];
562 fPoint2[0]=fTable->GetXCluster(1,ncl2);
563 fPoint2[1]=fTable->GetYCluster(1,ncl2);
565 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
566 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
567 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
568 trs,primaryVertex[2],pflag,fTable);
576 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
577 if(nn[nnp]!=0) fl+=1;
580 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
583 Int_t nct = trs->GetNumberOfClustersSA();
585 Int_t index = trs->GetClusterIndexSA(nct);
586 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
587 if(kl->IsUsed()==1) kl->Use();
593 AliESDtrack outtrack;
594 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
595 event->AddTrack(&outtrack);
596 Int_t nct = tr2->GetNumberOfClusters();
598 Int_t index = tr2->GetClusterIndex(nct);
599 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
601 if(kl->IsUsed()==0) kl->Use();
605 Int_t nct = trs->GetNumberOfClustersSA();
607 Int_t index = trs->GetClusterIndexSA(nct);
608 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
610 if(kl->IsUsed()==1) kl->Use();
614 }//end loop on clusters of layer2
621 Info("FindTracks","Number of found tracks: %d",event->GetNumberOfTracks());
628 //________________________________________________________________________
630 AliITStrackV2* AliITStrackerSA::FitTrack(AliITStrackSA* tr,Double_t *primaryVertex,Double_t *errorsprimvert){
631 //fit of the found track
634 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
635 for(Int_t i=0;i<fGeom->GetNlayers();i++){
636 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
638 AliITStrackV2* otrack2;
639 Int_t nclusters = tr->GetNumberOfClustersSA();
640 TObjArray** listlayer = new TObjArray*[fGeom->GetNlayers()];
641 for(Int_t i=0;i<fGeom->GetNlayers();i++){
642 listlayer[i] = new TObjArray(0,0);
652 Int_t * nnn = new Int_t[fGeom->GetNlayers()];
653 for(Int_t i=0;i<fGeom->GetNlayers();i++)nnn[i]=0;
655 for(Int_t ncl=0;ncl<nclusters;ncl++){
656 Int_t index = tr->GetClusterIndexSA(ncl);
657 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(index);
659 if(cl->IsUsed()==1) cl->Use();
660 Int_t lay = (index & 0xf0000000) >> 28;
661 if(lay==0) { listlayer[0]->AddLast(cl); clind0[nnn[0]]=index;nnn[0]++;}
662 if(lay==1) { listlayer[1]->AddLast(cl); clind1[nnn[1]]=index;nnn[1]++;}
663 if(lay==2) { listlayer[2]->AddLast(cl); clind2[nnn[2]]=index;nnn[2]++;}
664 if(lay==3) { listlayer[3]->AddLast(cl); clind3[nnn[3]]=index;nnn[3]++;}
665 if(lay==4) { listlayer[4]->AddLast(cl); clind4[nnn[4]]=index;nnn[4]++;}
666 if(lay==5) { listlayer[5]->AddLast(cl); clind5[nnn[5]]=index;nnn[5]++;}
671 Int_t * end = new Int_t[fGeom->GetNlayers()];
672 for(Int_t i=0;i<fGeom->GetNlayers();i++){
673 if(listlayer[i]->GetEntries()==0) end[i]=1;
674 else end[i]=listlayer[i]->GetEntries();
677 for(Int_t l1=0;l1<end[0];l1++){//loop on layer 1
678 AliITSclusterV2* cl0 = (AliITSclusterV2*)listlayer[0]->At(l1);
679 TVector3** recp = new TVector3*[3];
680 TVector3** errs = new TVector3*[3];
681 recp[0] = new TVector3(primaryVertex[0],primaryVertex[1],primaryVertex[2]);
682 errs[0] = new TVector3(errorsprimvert[0],errorsprimvert[1],errorsprimvert[2]);
683 Double_t x1,y1,z1,sx1,sy1,sz1;
684 Double_t x2,y2,z2,sx2,sy2,sz2;
685 AliITSclusterV2* p1=0;
686 AliITSclusterV2* p2=0;
687 Int_t index1=clind0[l1];
689 for(Int_t l2=0;l2<end[1];l2++){//loop on layer 2
690 AliITSclusterV2* cl1 = (AliITSclusterV2*)listlayer[1]->At(l2);
692 for(Int_t l3=0;l3<end[2];l3++){ //loop on layer 3
693 AliITSclusterV2* cl2 = (AliITSclusterV2*)listlayer[2]->At(l3);
695 if(cl0==0 && cl1!=0) {
696 p2 = cl2;index1=clind2[l3];
700 if(cl0!=0 && cl1==0){
702 p2=cl2;index2=clind2[l3];
704 if(cl0!=0 && cl1!=0){
708 Int_t lay1=(index1 & 0xf0000000) >> 28;
709 Int_t cln1=(index1 & 0x0fffffff) >> 00;
710 Int_t lay2=(index2 & 0xf0000000) >> 28;
711 Int_t cln2=(index2 & 0x0fffffff) >> 00;
712 x1 = fTable->GetXCluster(lay1,cln1);
713 x2 = fTable->GetXCluster(lay2,cln2);
714 y1 = fTable->GetYCluster(lay1,cln1);
715 y2 = fTable->GetYCluster(lay2,cln2);
716 z1 = fTable->GetZCluster(lay1,cln1);
717 z2 = fTable->GetZCluster(lay2,cln2);
718 sx1 = fTable->GetXClusterError(lay1,cln1);
719 sx2 = fTable->GetXClusterError(lay2,cln2);
720 sy1 = fTable->GetYClusterError(lay1,cln1);
721 sy2 = fTable->GetYClusterError(lay2,cln2);
722 sz1 = fTable->GetZClusterError(lay1,cln1);
723 sz2 = fTable->GetZClusterError(lay2,cln2);
724 Double_t phi1 = fTable->GetPhiCluster(lay1,cln1);
725 Int_t module1 = p1->GetDetectorIndex()+firstmod[0];
726 recp[1] = new TVector3(x1,y1,z1);
727 errs[1] = new TVector3(sx1,sy1,sz1);
728 recp[2] = new TVector3(x2,y2,z2);
729 errs[2] = new TVector3(sx2,sy2,sz2);
731 //fit on the Riemann sphere
732 Float_t seed1,seed2,seed3;
733 AliITSRiemannFit fit;
734 Int_t rf = fit.FitHelix(3,recp,errs,seed1,seed2,seed3); //this gives phi,tgl,curvature to start Kalman Filter
740 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
741 phi=seed1+1.5*TMath::Pi();
744 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
745 phi=seed1+0.5*TMath::Pi();
748 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
749 phi=seed1-0.5*TMath::Pi();
754 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
755 phi=seed1+0.5*TMath::Pi();
758 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
759 phi=seed1-0.5*TMath::Pi();
762 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
763 phi=seed1-1.5*TMath::Pi();
768 Int_t layer,ladder,detector;
769 fGeom->GetModuleId(module1,layer,ladder,detector);
770 Float_t yclu1 = p1->GetY();
771 Float_t zclu1 = p1->GetZ();
772 Double_t cv=Curvature(primaryVertex[0],primaryVertex[1],x1,y1,x2,y2);
774 for(Int_t l4=0;l4<end[3];l4++){ //loop on layer 4
775 AliITSclusterV2* cl3 = (AliITSclusterV2*)listlayer[3]->At(l4);
776 for(Int_t l5=0;l5<end[4];l5++){ //loop on layer 5
777 AliITSclusterV2* cl4 = (AliITSclusterV2*)listlayer[4]->At(l5);
778 for(Int_t l6=0;l6<end[5];l6++){ //loop on layer 6
779 AliITSclusterV2* cl5 = (AliITSclusterV2*)listlayer[5]->At(l6);
780 AliITStrackSA* trac = new AliITStrackSA(layer,ladder,detector,yclu1,zclu1,phi,tgl,cv,1);
782 if(cl5!=0) trac->AddClusterV2(5,(clind5[l6] & 0x0fffffff)>>0);
783 if(cl4!=0) trac->AddClusterV2(4,(clind4[l5] & 0x0fffffff)>>0);
784 if(cl3!=0) trac->AddClusterV2(3,(clind3[l4] & 0x0fffffff)>>0);
785 if(cl2!=0) trac->AddClusterV2(2,(clind2[l3] & 0x0fffffff)>>0);
786 if(cl1!=0) trac->AddClusterV2(1,(clind1[l2] & 0x0fffffff)>>0);
787 if(cl0!=0) trac->AddClusterV2(0,(clind0[l1] & 0x0fffffff)>>0);
789 //fit with Kalman filter using AliITStrackerV2::RefitAt()
791 AliITStrackV2* ot = new AliITStrackV2(*trac);
793 ot->ResetCovariance();
796 if(RefitAt(49.,ot,trac)){ //fit from layer 1 to layer 6
798 otrack2 = new AliITStrackV2(*ot);
799 otrack2->ResetCovariance();
800 otrack2->ResetClusters();
801 //fit from layer 6 to layer 1
802 if(RefitAt(3.7,otrack2,ot)) fListOfTracks->AddLast(otrack2);
812 for(Int_t i=1;i<3;i++){
826 Int_t dim=fListOfTracks->GetEntries();
828 for(Int_t i=0;i<fGeom->GetNlayers();i++){
835 AliITStrackV2* otrack =(AliITStrackV2*)FindTrackLowChiSquare(fListOfTracks,dim);
837 if(otrack==0) return 0;
838 Int_t * indexc = new Int_t[fGeom->GetNlayers()];
839 for(Int_t i=0;i<fGeom->GetNlayers();i++) indexc[i]=0;
840 for(Int_t nind=0;nind<otrack->GetNumberOfClusters();nind++){
841 indexc[nind] = otrack->GetClusterIndex(nind);
843 AliITSclusterV2* cl0 = (AliITSclusterV2*)GetCluster(indexc[0]);
844 AliITSclusterV2* cl1 = (AliITSclusterV2*)GetCluster(indexc[1]);
845 AliITSclusterV2* cl2 = (AliITSclusterV2*)GetCluster(indexc[2]);
846 AliITSclusterV2* cl3 = (AliITSclusterV2*)GetCluster(indexc[3]);
847 AliITSclusterV2* cl4 = (AliITSclusterV2*)GetCluster(indexc[4]);
848 Int_t labl[3]={-1,-1,-1};
849 if(otrack->GetNumberOfClusters()==fGeom->GetNlayers()){
850 AliITSclusterV2* cl5 = (AliITSclusterV2*)GetCluster(indexc[5]);
851 labl[0]=cl5->GetLabel(0);
852 labl[1]=cl5->GetLabel(1);
853 labl[2]=cl5->GetLabel(2);
856 if(otrack->GetNumberOfClusters()==(fGeom->GetNlayers()-1)){
861 Int_t numberofpoints;
862 if(fSixPoints) numberofpoints=6;
863 else numberofpoints=5;
864 Int_t label = Label(cl0->GetLabel(0),cl1->GetLabel(0),
865 cl2->GetLabel(0),cl3->GetLabel(0),
866 cl4->GetLabel(0),labl[0],
867 cl0->GetLabel(1),cl1->GetLabel(1),
868 cl2->GetLabel(1),cl3->GetLabel(1),
869 cl4->GetLabel(1),labl[1],
870 cl0->GetLabel(2),cl1->GetLabel(2),
871 cl2->GetLabel(2),cl3->GetLabel(2),
872 cl4->GetLabel(2),labl[2],numberofpoints);
874 otrack->SetLabel(label);
875 for(Int_t i=0;i<fGeom->GetNlayers();i++){
884 //_______________________________________________________________________
885 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,TTree* treev2){
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)
892 if(fVert)delete fVert;
893 fVert = fVertexer->FindVertexForCurrentEvent(evnum);
896 gAlice->GetEvent(evnum);
898 Fatal("FindTracks","Vertex is missing\n");
902 Double_t primaryVertex[3];
903 fVert->GetXYZ(primaryVertex);
906 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
907 fTable->FillArray(fITSclusters);
908 fTable->FillArrayCoorAngles();
911 TBranch* bra = (TBranch*)treev2->GetBranch("tracks");
912 if(!bra) Warning("UseFoundTracksV2","No branch for track tree");
913 AliITStrackV2* ttrrt = new AliITStrackV2;
914 bra->SetAddress(&ttrrt);
916 for(Int_t nj=0;nj<treev2->GetEntries();nj++){
917 treev2->GetEvent(nj);
918 Int_t ncl = ttrrt->GetNumberOfClusters();
919 for(Int_t k=0;k<ncl;k++){
920 Int_t index = ttrrt->GetClusterIndex(k);
921 AliITSclusterV2* clui = (AliITSclusterV2*)GetCluster(index);
922 if(clui->IsUsed()==0) clui->Use();
929 //_______________________________________________________________________
930 void AliITStrackerSA::UseFoundTracksV2(AliESD *event){
931 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
935 Double_t primaryVertex[3];
936 event->GetVertex()->GetXYZ(primaryVertex);
939 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
940 fTable->FillArray(fITSclusters);
941 fTable->FillArrayCoorAngles();
944 Int_t ntracks = event->GetNumberOfTracks();
946 AliESDtrack *esd=event->GetTrack(ntracks);
947 if ((esd->GetStatus()&
948 AliESDtrack::kITSin|AliESDtrack::kTPCin)==0) continue;
950 Int_t ncl = esd->GetITSclusters(idx);
951 for(Int_t clu=0; clu<ncl; clu++){
952 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(idx[clu]);
953 if(cl->IsUsed()==0) cl->Use();
957 Info("UseFoundTracksV2","Clusters of tracks prolonged from TPC deleted");
962 //_______________________________________________________
963 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Int_t pflag,AliITSclusterTable* table){
964 //function used to to find the clusters associated to the track
966 AliITSlayer &lay = fgLayers[layer];
967 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
968 for(Int_t i=0;i<fGeom->GetNlayers();i++){
969 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
973 Float_t cx1,cx2,cy1,cy2;
974 FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
975 Int_t fun = FindIntersection(fCoef1,fCoef2,fCoef3,-(lay.GetR()*lay.GetR()),cx1,cy1,cx2,cy2);
980 Double_t fi1 =TMath::ATan2(cy1,cx1);
981 Double_t fi2 =TMath::ATan2(cy2,cx2);
982 fPhiEstimate = ChoosePoint(fi1,fi2,fPhic);
985 Double_t zed = TMath::Tan(fLambdac)*lay.GetR()+zvertex;
986 Double_t zed1 = TMath::Tan(fLambdac+lambdawindow)*lay.GetR()+zvertex;
987 Double_t zed2 = TMath::Tan(fLambdac-lambdawindow)*lay.GetR()+zvertex;
989 Double_t fi = fPhiEstimate;
990 Int_t nmod = lay.FindDetectorIndex(fi,zed);
995 nmod += firstmod[layer];
997 Int_t nm[8]={0,0,0,0,0,0,0,0};
998 nm[0] = lay.FindDetectorIndex(fi+phiwindow,zed);
999 nm[1] = lay.FindDetectorIndex(fi-phiwindow,zed);
1000 nm[2] = lay.FindDetectorIndex(fi,zed1);
1001 nm[3] = lay.FindDetectorIndex(fi,zed2);
1002 nm[4] = lay.FindDetectorIndex(fi+phiwindow,zed1);
1003 nm[5] = lay.FindDetectorIndex(fi-phiwindow,zed1);
1004 nm[6] = lay.FindDetectorIndex(fi+phiwindow,zed2);
1005 nm[7] = lay.FindDetectorIndex(fi-phiwindow,zed2);
1009 TArrayI* array =(TArrayI*)table->GetListOfClusters(nmod);
1010 TArrayI* listc = new TArrayI(array->GetSize());
1011 for(Int_t i=0;i<array->GetSize();i++){
1012 Int_t in=(Int_t)array->At(i);
1013 listc->AddAt(in,nn);
1020 for(Int_t h=k+1;h<8;h++){
1033 for(Int_t ii=0;ii<8;ii++){
1034 if(nm[ii]!=value && nm[ii]!=nmod && nm[ii]>=0){
1035 TArrayI* ar =(TArrayI*)table->GetListOfClusters(nm[ii]+firstmod[layer]);
1036 listc->Set(listc->GetSize()+ar->GetSize());
1037 for(Int_t j=0;j<ar->GetSize();j++){
1038 Int_t in=(Int_t)ar->At(j);
1039 listc->AddAt(in,nn);
1047 for(Int_t i=0;i<listc->GetSize();i++){
1048 Int_t index = (Int_t)listc->At(i);
1049 AliITSclusterV2* cllay = lay.GetCluster(index);
1050 if(cllay==0) continue;
1051 if(cllay->IsUsed()==1) continue;
1052 Double_t phi = fTable->GetPhiCluster(layer,index);
1053 Double_t lambda= fTable->GetLambdaCluster(layer,index);
1055 if(TMath::Abs(fLambdac-lambda)<lambdawindow &&
1056 TMath::Abs(fPhiEstimate-phi)<phiwindow){
1059 if(trs->GetNumberOfClustersSA()==15){
1064 trs->AddClusterSA(layer,index);
1067 fPointc[0]=fTable->GetXCluster(layer,index);
1068 fPointc[1]=fTable->GetYCluster(layer,index);
1080 //________________________________________________________________
1081 void AliITStrackerSA::UpdatePoints(){
1082 //update of points for the estimation of the curvature
1084 fPoint1[0]=fPoint2[0];
1085 fPoint1[1]=fPoint2[1];
1086 fPoint2[0]=fPoint3[0];
1087 fPoint2[1]=fPoint3[1];
1088 fPoint3[0]=fPointc[0];
1089 fPoint3[1]=fPointc[1];
1096 //___________________________________________________________________
1097 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){
1099 //given (x,y) of three recpoints (in global coordinates)
1100 //returns the parameters a,b,c of circonference x*x + y*y +a*x + b*y +c
1102 Float_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1103 if(den==0) return 0;
1104 a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1105 b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1106 c = -x1*x1-y1*y1-a*x1-b*y1;
1109 //__________________________________________________________________________
1110 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){
1112 //Finds the intersection between the circonference of the track and the circonference centered in (0,0) represented by one layer
1113 //c2 is -rlayer*rlayer
1117 Float_t aA = (b1*b1)/(a1*a1)+1;
1118 Float_t bB = (-2*m*b1/(a1*a1));
1119 Float_t cC = c2+(m*m)/(a1*a1);
1120 if((bB*bB-4*aA*cC)<0) return 0;
1122 y1 = (-bB+TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA);
1123 y2 = (-bB-TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA);
1124 x1 = (c2-c1-b1*y1)/a1;
1125 x2 = (c2-c1-b1*y2)/a1;
1129 //____________________________________________________________________
1130 Double_t AliITStrackerSA::Curvature(Double_t x1,Double_t y1,Double_t
1131 x2,Double_t y2,Double_t x3,Double_t y3){
1133 //calculates the curvature of track
1134 Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1135 if(den==0) return 0;
1136 Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1137 Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1138 Double_t c = -x1*x1-y1*y1-a*x1-b*y1;
1141 if((a*a+b*b-4*c)<0) return 0;
1142 Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.;
1143 if(rad==0) return 0;
1145 if((x1>0 && y1>0 && x1<xc)) rad*=-1;
1146 if((x1<0 && y1>0 && x1<xc)) rad*=-1;
1147 // if((x1<0 && y1<0 && x1<xc)) rad*=-1;
1148 // if((x1>0 && y1<0 && x1<xc)) rad*=-1;
1153 //____________________________________________________________________
1154 Double_t AliITStrackerSA::ChoosePoint(Double_t p1, Double_t p2, Double_t pp){
1156 //Returns the point closest to pp
1158 Double_t diff1 = p1-pp;
1159 Double_t diff2 = p2-pp;
1161 if(TMath::Abs(diff1)<TMath::Abs(diff2)) fPhiEstimate=p1;
1162 else fPhiEstimate=p2;
1163 return fPhiEstimate;
1168 //_________________________________________________________________
1169 AliITStrackV2* AliITStrackerSA::FindTrackLowChiSquare(TObjArray* tracklist, Int_t dim) const {
1170 // returns track with lowes chi square
1172 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(0);
1175 if(dim==0) return 0;
1176 Double_t * chi2 = new Double_t[dim];
1177 Int_t * index = new Int_t[dim];
1178 for(Int_t i=0;i<dim;i++){
1179 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(i);
1180 chi2[i]=trk->GetChi2();
1184 Int_t w=0;Double_t value;
1187 for(Int_t j=w+1;j<dim;j++){
1188 if(chi2[w]<chi2[j]){
1200 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(index[dim-1]);
1207 //__________________________________________________________
1208 Int_t AliITStrackerSA::FindLabel(Int_t l1, Int_t l2, Int_t l3, Int_t l4, Int_t l5, Int_t l6){
1210 //function used to determine the track label
1212 Int_t lb[6] = {l1,l2,l3,l4,l5,l6};
1213 Int_t aa[6]={1,1,1,1,1,1};
1216 Int_t k=0;Int_t w=0;Int_t num=6;
1217 if(lb[5]==-1) num=5;
1221 for(Int_t i=k+1;i<num;i++){
1223 if(lb[k]==lb[i] && aa[k]!=0){
1234 for(Int_t j=0;j<6;j++){
1246 if(num==6) return lb[5];
1250 //_____________________________________________________________________________
1251 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,
1252 Int_t gl12, Int_t gl13, Int_t gl14,Int_t gl15, Int_t gl16, Int_t gl17, Int_t gl18, Int_t numberofpoints){
1255 //function used to assign label to the found track. If track is fake, the label is negative
1257 Int_t lb0[6] = {gl1,gl2,gl3,gl4,gl5,gl6};
1258 Int_t lb1[6] = {gl7,gl8,gl9,gl10,gl11,gl12};
1259 Int_t lb2[6] = {gl13,gl14,gl15,gl16,gl17,gl18};
1260 Int_t ll=FindLabel(lb0[0],lb0[1],lb0[2],lb0[3],lb0[4],lb0[5]);
1261 Int_t lflag=0;Int_t num=6;
1262 if(lb0[5]==-1 && lb1[5]==-1 && lb2[5]==-1) num=5;
1264 for(Int_t i=0;i<num;i++){
1265 if(lb0[i]==ll || lb1[i]==ll || lb2[i]==ll) lflag+=1;
1268 if(lflag>=numberofpoints) return ll;
1274 //_____________________________________________________________________________
1275 void AliITStrackerSA::SetWindowSizes(Int_t n, Double_t *phi, Double_t *lam){
1276 // Set sizes of the phi and lambda windows used for track finding
1278 if(phi){ // user defined values
1279 fPhiWin = new Double_t[fNloop];
1280 fLambdaWin = new Double_t[fNloop];
1281 for(Int_t k=0;k<fNloop;k++){
1283 fLambdaWin[k]=lam[k];
1286 else { // default values
1288 Double_t phid[32] = {0.002,0.003,0.004,0.0045,0.0047,
1289 0.005,0.0053,0.0055,
1290 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1291 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1292 0.01,0.0105,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1294 Double_t lambdad[32] = {0.002,0.003,0.004,0.0045,0.0047,
1295 0.005,0.0053,0.0055,
1296 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1297 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1298 0.01,0.015,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1305 fPhiWin = new Double_t[fNloop];
1306 fLambdaWin = new Double_t[fNloop];
1308 for(Int_t k=0;k<fNloop;k++){
1310 fLambdaWin[k]=lambdad[k];