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 ////////////////////////////////////////////////////
28 #include <TObjArray.h>
31 #include "AliITSclusterTable.h"
32 #include "AliITSclusterV2.h"
33 #include "AliITSgeom.h"
34 #include "AliITSRiemannFit.h"
35 #include "AliITStrackerSA.h"
36 #include "AliITStrackSA.h"
37 #include "AliITSVertexer.h"
38 #include "AliESDVertex.h"
40 #include "AliESDtrack.h"
42 ClassImp(AliITStrackerSA)
44 //____________________________________________________________________________
45 AliITStrackerSA::AliITStrackerSA():AliITStrackerV2(){
46 // Default constructor
49 //____________________________________________________________________________
50 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom):AliITStrackerV2(geom)
52 // Standard constructor (Vertex is known and passed to this obj.)
58 //____________________________________________________________________________
59 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliESDVertex *vert):AliITStrackerV2(geom)
61 // Standard constructor (Vertex is known and passed to this obj.)
67 //______________________________________________________________________
68 AliITStrackerSA::AliITStrackerSA(const AliITStrackerSA &trkr) :
69 AliITStrackerV2(trkr) {
71 // Copies are not allowed. The method is protected to avoid misuse.
72 Error("AliITStrackerSA","Copy constructor not allowed\n");
75 //______________________________________________________________________
76 AliITStrackerSA& AliITStrackerSA::operator=(const
77 AliITStrackerSA& /* trkr */){
78 // Assignment operator
79 // Assignment is not allowed. The method is protected to avoid misuse.
80 Error("= operator","Assignment operator not allowed\n");
84 //____________________________________________________________________________
85 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertexer *vertexer):AliITStrackerV2(geom)
87 // Standard constructor (Vertex is unknown - vertexer is passed to this obj)
93 //____________________________________________________________________________
94 AliITStrackerSA::AliITStrackerSA(AliITStrackerSA& tracker):AliITStrackerV2(){
96 fPhiEstimate = tracker.fPhiEstimate;
97 for(Int_t i=0;i<2;i++){
98 fPoint1[i]=tracker.fPoint1[i];
99 fPoint2[i]=tracker.fPoint2[i];
100 fPoint3[i]=tracker.fPoint3[i];
101 fPointc[i]=tracker.fPointc[i];
103 fLambdac = tracker.fLambdac;
104 fPhic = tracker.fPhic;
105 fCoef1 = tracker.fCoef1;
106 fCoef2 = tracker.fCoef2;
107 fCoef3 = tracker.fCoef3;
108 fNloop = tracker.fNloop;
109 fPhiWin = tracker.fPhiWin;
110 fLambdaWin = tracker.fLambdaWin;
111 if(tracker.fVertexer && tracker.fVert){
112 fVert = new AliESDVertex(*tracker.fVert);
115 fVert = tracker.fVert;
117 fVertexer = tracker.fVertexer;
118 fGeom = tracker.fGeom;
119 fTable = tracker.fTable;
120 fListOfTracks = tracker.fListOfTracks;
123 //____________________________________________________________________________
124 AliITStrackerSA::~AliITStrackerSA(){
126 // if fVertexer is not null, the AliESDVertex obj. is owned by this class
127 // and is deleted here
129 if(fVert)delete fVert;
134 if(fPhiWin)delete []fPhiWin;
135 if(fLambdaWin)delete []fLambdaWin;
137 fListOfTracks->Delete();
140 //____________________________________________________________________________
141 void AliITStrackerSA::Init(){
142 // Reset all data members
144 for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
159 fListOfTracks=new TObjArray(0,0);
161 //_______________________________________________________________________
162 void AliITStrackerSA::ResetForFinding(){
163 // Reset data members used in all loops during track finding
165 for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
173 fListOfTracks->Delete();
175 //____________________________________________________________________________
176 void AliITStrackerSA::FindTracks(TTree *out,Int_t evnumber){
178 /**************************************************************************
179 * This function finds primary tracks.
182 * Example: to execute function with only the ITS (no combined tracking *
183 * with TPC+ITS) and requiring 5/6 points to define a good track *
184 * call SetSixPoinbts(kFALSE) in advance and then *
185 * use: FindTracks(treein,treeout,evnumber) *
186 * to execute combined tracking, before using FindTracks, use *
188 *************************************************************************/
191 Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
196 if(fVert)delete fVert;
197 fVert = fVertexer->FindVertexForCurrentEvent(evnumber);
200 gAlice->GetEvent(evnumber);
202 Fatal("FindTracks","Vertex is missing\n");
206 Double_t primaryVertex[3];
207 Double_t errorsprimvert[3];
208 fVert->GetXYZ(primaryVertex);
209 fVert->GetSigmaXYZ(errorsprimvert);
210 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
211 Warning("FindTracks","Set errors on vertex positions x and y at 0.0001");
212 errorsprimvert[0]=0.0001;
213 errorsprimvert[1]=0.0001;
215 fVert->PrintStatus();
218 //Fill array with cluster indices for each module
220 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
221 fTable->FillArray(fITSclusters);
222 fTable->FillArrayCoorAngles();
226 //Fill tree for found tracks
227 AliITStrackV2* outrack=0;
228 TBranch* branch=out->Branch("tracks","AliITStrackV2",&outrack,32000,0);
229 if (!branch) out->Branch("tracks","AliITStrackV2",&outrack,32000,3);
230 else branch->SetAddress(&outrack);
233 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
234 for(Int_t i=0;i<fGeom->GetNlayers();i++){
235 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
237 // firstmod [i] number of the first module in the ITS layer i.
239 AliITSlayer &layer=fgLayers[0]; // first layer
241 Int_t dim=layer.GetNumberOfClusters();
242 //loop on the different windows
243 for(Int_t nloop=0;nloop<fNloop;nloop++){
244 for(Int_t ncl=0;ncl<dim;ncl++){
245 //loop starting from layer 0
248 AliITSclusterV2* cl = layer.GetCluster(ncl);
249 if(cl->IsUsed()==1) continue;
250 fPhic = fTable->GetPhiCluster(0,ncl);
251 fLambdac = fTable->GetLambdaCluster(0,ncl);
252 fPhiEstimate = fPhic;
253 AliITStrackSA* trs = new AliITStrackSA();
254 fPoint1[0]=primaryVertex[0];
255 fPoint1[1]=primaryVertex[1];
256 fPoint2[0]=fTable->GetXCluster(0,ncl);
257 fPoint2[1]=fTable->GetYCluster(0,ncl);
259 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
260 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
261 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
262 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
266 fPoint3[0] = fPointc[0];
267 fPoint3[1] = fPointc[1];
269 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
270 if(nn[1]==0 && nn[2]==0) pflag=0;
271 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
272 if(nn[2]!=0 && nn[1]==0){
274 fPoint3[0]=fPointc[0];
275 fPoint3[1]=fPointc[1];
278 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
280 if(nn[3]!=0) UpdatePoints();
281 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
283 if(nn[4]!=0) UpdatePoints();
284 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
288 Int_t numberofpoints;
289 if(fSixPoints) numberofpoints=6; //check of the candidate track
290 else numberofpoints=5; //if track is good (with the required number
291 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
292 if(nn[nnp]!=0) layOK+=1;
294 if(layOK>=numberofpoints){
295 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
297 Int_t nct = trs->GetNumberOfClustersSA();
299 Int_t index = trs->GetClusterIndexSA(nct);
300 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
301 if(kl->IsUsed()==1) kl->Use();
308 Int_t nct = tr2->GetNumberOfClusters();
311 Int_t index = tr2->GetClusterIndex(nct);
312 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
313 if(kl->IsUsed()==0) kl->Use();
318 Int_t nct = trs->GetNumberOfClustersSA();
320 Int_t index = trs->GetClusterIndexSA(nct);
321 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
322 if(kl->IsUsed()==1) kl->Use();
327 }//end loop on clusters of layer1
331 //if 5/6 points are required, second loop starting
332 //from second layer, to find tracks with point of
336 // counter for clusters on each layer
337 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
338 for(Int_t nloop=0;nloop<fNloop;nloop++){
339 AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
340 Int_t ncl2=layer2.GetNumberOfClusters();
341 while(ncl2--){ //loop starting from layer 2
344 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
345 if(cl->IsUsed()==1) continue;
346 fPhic = fTable->GetPhiCluster(1,ncl2);
347 fLambdac = fTable->GetLambdaCluster(1,ncl2);
348 fPhiEstimate = fPhic;
349 AliITStrackSA* trs = new AliITStrackSA();
351 fPoint1[0]=primaryVertex[0];
352 fPoint1[1]=primaryVertex[1];
353 fPoint2[0]=fTable->GetXCluster(1,ncl2);;
354 fPoint2[1]=fTable->GetYCluster(1,ncl2);;
356 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
357 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
358 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
359 trs,primaryVertex[2],pflag,fTable);
367 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
368 if(nn[nnp]!=0) fl+=1;
371 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
373 Int_t nct = trs->GetNumberOfClustersSA();
375 Int_t index = trs->GetClusterIndexSA(nct);
376 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
377 if(kl->IsUsed()==1) kl->Use();
383 Int_t nct = tr2->GetNumberOfClusters();
385 Int_t index = tr2->GetClusterIndex(nct);
386 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
388 if(kl->IsUsed()==0) kl->Use();
392 Int_t nct = trs->GetNumberOfClustersSA();
394 Int_t index = trs->GetClusterIndexSA(nct);
395 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
397 if(kl->IsUsed()==1) kl->Use();
401 }//end loop on clusters of layer2
404 } // if(!fSixPoints....
412 //______________________________________________________________________
413 Int_t AliITStrackerSA::FindTracks(AliESD* event){
415 // Track finder using the ESD object
419 Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
424 Double_t errorsprimvert[3];
425 Double_t primaryVertex[3];
426 event->GetVertex()->GetXYZ(primaryVertex);
427 event->GetVertex()->GetSigmaXYZ(errorsprimvert);
429 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
430 // Warning("FindTracks","Set errors on vertex positions x and y at 0.005");
431 errorsprimvert[0]=0.005;
432 errorsprimvert[1]=0.005;
435 //Fill array with cluster indices for each module
437 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
438 fTable->FillArray(fITSclusters);
439 fTable->FillArrayCoorAngles();
442 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
443 for(Int_t i=0;i<fGeom->GetNlayers();i++){
444 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
446 // firstmod [i] number of the first module in the ITS layer i.
449 AliITSlayer &layer=fgLayers[0]; // first layer
452 Int_t dim=layer.GetNumberOfClusters();
453 //loop on the different windows
454 for(Int_t nloop=0;nloop<fNloop;nloop++){
455 for(Int_t ncl=0;ncl<dim;ncl++){ //loop starting from layer 0
458 AliITSclusterV2* cl = (AliITSclusterV2*)layer.GetCluster(ncl);
460 if(cl->IsUsed()==1) continue;
461 fPhic = fTable->GetPhiCluster(0,ncl);
462 fLambdac = fTable->GetLambdaCluster(0,ncl);
463 fPhiEstimate = fPhic;
464 AliITStrackSA* trs = new AliITStrackSA();
465 fPoint1[0]=primaryVertex[0];
466 fPoint1[1]=primaryVertex[1];
467 fPoint2[0]=fTable->GetXCluster(0,ncl);
468 fPoint2[1]=fTable->GetYCluster(0,ncl);
469 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
470 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
471 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
473 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
476 fPoint3[0] = fPointc[0];
477 fPoint3[1] = fPointc[1];
479 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
480 if(nn[1]==0 && nn[2]==0) pflag=0;
481 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
482 if(nn[2]!=0 && nn[1]==0){
484 fPoint3[0]=fPointc[0];
485 fPoint3[1]=fPointc[1];
488 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
490 if(nn[3]!=0) UpdatePoints();
491 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
493 if(nn[4]!=0) UpdatePoints();
494 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
498 Int_t numberofpoints;
499 if(fSixPoints) numberofpoints=6; //check of the candidate track
500 else numberofpoints=5; //if track is good (with the required number
501 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
502 if(nn[nnp]!=0) layOK+=1;
504 if(layOK>=numberofpoints){
505 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
507 Int_t nct = trs->GetNumberOfClustersSA();
509 Int_t index = trs->GetClusterIndexSA(nct);
510 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
511 if(kl->IsUsed()==1) kl->Use();
516 AliESDtrack outtrack;
517 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
518 event->AddTrack(&outtrack);
520 Int_t nct = tr2->GetNumberOfClusters();
522 Int_t index = tr2->GetClusterIndex(nct);
523 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
524 if(kl->IsUsed()==0) kl->Use();
529 Int_t nct = trs->GetNumberOfClustersSA();
531 Int_t index = trs->GetClusterIndexSA(nct);
532 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
533 if(kl->IsUsed()==1) kl->Use();
539 }//end loop on clusters of layer1
543 //if 5/6 points are required, second loop starting
544 //from second layer, to find tracks with point of
548 // counter for clusters on each layer
549 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
550 for(Int_t nloop=0;nloop<fNloop;nloop++){
551 AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
552 Int_t ncl2=layer2.GetNumberOfClusters();
553 while(ncl2--){ //loop starting from layer 2
556 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
557 if(cl->IsUsed()==1) continue;
558 fPhic = fTable->GetPhiCluster(1,ncl2);
559 fLambdac = fTable->GetLambdaCluster(1,ncl2);
560 fPhiEstimate = fPhic;
561 AliITStrackSA* trs = new AliITStrackSA();
562 fPoint1[0]=primaryVertex[0];
563 fPoint1[1]=primaryVertex[1];
564 fPoint2[0]=fTable->GetXCluster(1,ncl2);
565 fPoint2[1]=fTable->GetYCluster(1,ncl2);
567 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
568 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
569 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
570 trs,primaryVertex[2],pflag,fTable);
578 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
579 if(nn[nnp]!=0) fl+=1;
582 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
585 Int_t nct = trs->GetNumberOfClustersSA();
587 Int_t index = trs->GetClusterIndexSA(nct);
588 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
589 if(kl->IsUsed()==1) kl->Use();
595 AliESDtrack outtrack;
596 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
597 event->AddTrack(&outtrack);
598 Int_t nct = tr2->GetNumberOfClusters();
600 Int_t index = tr2->GetClusterIndex(nct);
601 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
603 if(kl->IsUsed()==0) kl->Use();
607 Int_t nct = trs->GetNumberOfClustersSA();
609 Int_t index = trs->GetClusterIndexSA(nct);
610 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
612 if(kl->IsUsed()==1) kl->Use();
616 }//end loop on clusters of layer2
624 Info("FindTracks","Number of found tracks: %d",event->GetNumberOfTracks());
631 //________________________________________________________________________
633 AliITStrackV2* AliITStrackerSA::FitTrack(AliITStrackSA* tr,Double_t *primaryVertex,Double_t *errorsprimvert){
634 //fit of the found track
637 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
638 for(Int_t i=0;i<fGeom->GetNlayers();i++){
639 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
641 AliITStrackV2* otrack2;
642 Int_t nclusters = tr->GetNumberOfClustersSA();
643 TObjArray** listlayer = new TObjArray*[fGeom->GetNlayers()];
644 for(Int_t i=0;i<fGeom->GetNlayers();i++){
645 listlayer[i] = new TObjArray(0,0);
655 Int_t * nnn = new Int_t[fGeom->GetNlayers()];
656 for(Int_t i=0;i<fGeom->GetNlayers();i++)nnn[i]=0;
658 for(Int_t ncl=0;ncl<nclusters;ncl++){
659 Int_t index = tr->GetClusterIndexSA(ncl);
660 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(index);
662 if(cl->IsUsed()==1) cl->Use();
663 Int_t lay = (index & 0xf0000000) >> 28;
664 if(lay==0) { listlayer[0]->AddLast(cl); clind0[nnn[0]]=index;nnn[0]++;}
665 if(lay==1) { listlayer[1]->AddLast(cl); clind1[nnn[1]]=index;nnn[1]++;}
666 if(lay==2) { listlayer[2]->AddLast(cl); clind2[nnn[2]]=index;nnn[2]++;}
667 if(lay==3) { listlayer[3]->AddLast(cl); clind3[nnn[3]]=index;nnn[3]++;}
668 if(lay==4) { listlayer[4]->AddLast(cl); clind4[nnn[4]]=index;nnn[4]++;}
669 if(lay==5) { listlayer[5]->AddLast(cl); clind5[nnn[5]]=index;nnn[5]++;}
674 Int_t * end = new Int_t[fGeom->GetNlayers()];
675 for(Int_t i=0;i<fGeom->GetNlayers();i++){
676 if(listlayer[i]->GetEntries()==0) end[i]=1;
677 else end[i]=listlayer[i]->GetEntries();
680 for(Int_t l1=0;l1<end[0];l1++){//loop on layer 1
681 AliITSclusterV2* cl0 = (AliITSclusterV2*)listlayer[0]->At(l1);
682 TVector3** recp = new TVector3*[3];
683 TVector3** errs = new TVector3*[3];
684 recp[0] = new TVector3(primaryVertex[0],primaryVertex[1],primaryVertex[2]);
685 errs[0] = new TVector3(errorsprimvert[0],errorsprimvert[1],errorsprimvert[2]);
686 Double_t x1,y1,z1,sx1,sy1,sz1;
687 Double_t x2,y2,z2,sx2,sy2,sz2;
688 AliITSclusterV2* p1=0;
689 AliITSclusterV2* p2=0;
690 Int_t index1=clind0[l1];
692 for(Int_t l2=0;l2<end[1];l2++){//loop on layer 2
693 AliITSclusterV2* cl1 = (AliITSclusterV2*)listlayer[1]->At(l2);
695 for(Int_t l3=0;l3<end[2];l3++){ //loop on layer 3
696 AliITSclusterV2* cl2 = (AliITSclusterV2*)listlayer[2]->At(l3);
698 if(cl0==0 && cl1!=0) {
699 p2 = cl2;index1=clind2[l3];
703 if(cl0!=0 && cl1==0){
705 p2=cl2;index2=clind2[l3];
707 if(cl0!=0 && cl1!=0){
711 Int_t lay1=(index1 & 0xf0000000) >> 28;
712 Int_t cln1=(index1 & 0x0fffffff) >> 00;
713 Int_t lay2=(index2 & 0xf0000000) >> 28;
714 Int_t cln2=(index2 & 0x0fffffff) >> 00;
715 x1 = fTable->GetXCluster(lay1,cln1);
716 x2 = fTable->GetXCluster(lay2,cln2);
717 y1 = fTable->GetYCluster(lay1,cln1);
718 y2 = fTable->GetYCluster(lay2,cln2);
719 z1 = fTable->GetZCluster(lay1,cln1);
720 z2 = fTable->GetZCluster(lay2,cln2);
721 sx1 = fTable->GetXClusterError(lay1,cln1);
722 sx2 = fTable->GetXClusterError(lay2,cln2);
723 sy1 = fTable->GetYClusterError(lay1,cln1);
724 sy2 = fTable->GetYClusterError(lay2,cln2);
725 sz1 = fTable->GetZClusterError(lay1,cln1);
726 sz2 = fTable->GetZClusterError(lay2,cln2);
727 Double_t phi1 = fTable->GetPhiCluster(lay1,cln1);
728 Int_t module1 = p1->GetDetectorIndex()+firstmod[0];
729 recp[1] = new TVector3(x1,y1,z1);
730 errs[1] = new TVector3(sx1,sy1,sz1);
731 recp[2] = new TVector3(x2,y2,z2);
732 errs[2] = new TVector3(sx2,sy2,sz2);
734 //fit on the Riemann sphere
735 Float_t seed1,seed2,seed3;
736 AliITSRiemannFit fit;
737 Int_t rf = fit.FitHelix(3,recp,errs,seed1,seed2,seed3); //this gives phi,tgl,curvature to start Kalman Filter
743 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
744 phi=seed1+1.5*TMath::Pi();
747 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
748 phi=seed1+0.5*TMath::Pi();
751 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
752 phi=seed1-0.5*TMath::Pi();
757 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
758 phi=seed1+0.5*TMath::Pi();
761 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
762 phi=seed1-0.5*TMath::Pi();
765 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
766 phi=seed1-1.5*TMath::Pi();
771 Int_t layer,ladder,detector;
772 fGeom->GetModuleId(module1,layer,ladder,detector);
773 Float_t yclu1 = p1->GetY();
774 Float_t zclu1 = p1->GetZ();
775 Double_t cv=Curvature(primaryVertex[0],primaryVertex[1],x1,y1,x2,y2);
777 for(Int_t l4=0;l4<end[3];l4++){ //loop on layer 4
778 AliITSclusterV2* cl3 = (AliITSclusterV2*)listlayer[3]->At(l4);
779 for(Int_t l5=0;l5<end[4];l5++){ //loop on layer 5
780 AliITSclusterV2* cl4 = (AliITSclusterV2*)listlayer[4]->At(l5);
781 for(Int_t l6=0;l6<end[5];l6++){ //loop on layer 6
782 AliITSclusterV2* cl5 = (AliITSclusterV2*)listlayer[5]->At(l6);
783 AliITStrackSA* trac = new AliITStrackSA(layer,ladder,detector,yclu1,zclu1,phi,tgl,cv,1);
785 if(cl5!=0) trac->AddClusterV2(5,(clind5[l6] & 0x0fffffff)>>0);
786 if(cl4!=0) trac->AddClusterV2(4,(clind4[l5] & 0x0fffffff)>>0);
787 if(cl3!=0) trac->AddClusterV2(3,(clind3[l4] & 0x0fffffff)>>0);
788 if(cl2!=0) trac->AddClusterV2(2,(clind2[l3] & 0x0fffffff)>>0);
789 if(cl1!=0) trac->AddClusterV2(1,(clind1[l2] & 0x0fffffff)>>0);
790 if(cl0!=0) trac->AddClusterV2(0,(clind0[l1] & 0x0fffffff)>>0);
792 //fit with Kalman filter using AliITStrackerV2::RefitAt()
794 AliITStrackV2* ot = new AliITStrackV2(*trac);
796 ot->ResetCovariance();
799 if(RefitAt(49.,ot,trac)){ //fit from layer 1 to layer 6
801 otrack2 = new AliITStrackV2(*ot);
802 otrack2->ResetCovariance();
803 otrack2->ResetClusters();
804 //fit from layer 6 to layer 1
805 if(RefitAt(3.7,otrack2,ot)) fListOfTracks->AddLast(otrack2);
815 for(Int_t i=1;i<3;i++){
829 Int_t dim=fListOfTracks->GetEntries();
831 for(Int_t i=0;i<fGeom->GetNlayers();i++){
838 AliITStrackV2* otrack =(AliITStrackV2*)FindTrackLowChiSquare(fListOfTracks,dim);
840 if(otrack==0) return 0;
841 Int_t * indexc = new Int_t[fGeom->GetNlayers()];
842 for(Int_t i=0;i<fGeom->GetNlayers();i++) indexc[i]=0;
843 for(Int_t nind=0;nind<otrack->GetNumberOfClusters();nind++){
844 indexc[nind] = otrack->GetClusterIndex(nind);
846 AliITSclusterV2* cl0 = (AliITSclusterV2*)GetCluster(indexc[0]);
847 AliITSclusterV2* cl1 = (AliITSclusterV2*)GetCluster(indexc[1]);
848 AliITSclusterV2* cl2 = (AliITSclusterV2*)GetCluster(indexc[2]);
849 AliITSclusterV2* cl3 = (AliITSclusterV2*)GetCluster(indexc[3]);
850 AliITSclusterV2* cl4 = (AliITSclusterV2*)GetCluster(indexc[4]);
851 Int_t labl[3]={-1,-1,-1};
852 if(otrack->GetNumberOfClusters()==fGeom->GetNlayers()){
853 AliITSclusterV2* cl5 = (AliITSclusterV2*)GetCluster(indexc[5]);
854 labl[0]=cl5->GetLabel(0);
855 labl[1]=cl5->GetLabel(1);
856 labl[2]=cl5->GetLabel(2);
859 if(otrack->GetNumberOfClusters()==(fGeom->GetNlayers()-1)){
864 Int_t numberofpoints;
865 if(fSixPoints) numberofpoints=6;
866 else numberofpoints=5;
867 Int_t label = Label(cl0->GetLabel(0),cl1->GetLabel(0),
868 cl2->GetLabel(0),cl3->GetLabel(0),
869 cl4->GetLabel(0),labl[0],
870 cl0->GetLabel(1),cl1->GetLabel(1),
871 cl2->GetLabel(1),cl3->GetLabel(1),
872 cl4->GetLabel(1),labl[1],
873 cl0->GetLabel(2),cl1->GetLabel(2),
874 cl2->GetLabel(2),cl3->GetLabel(2),
875 cl4->GetLabel(2),labl[2],numberofpoints);
877 otrack->SetLabel(label);
878 for(Int_t i=0;i<fGeom->GetNlayers();i++){
887 //_______________________________________________________________________
888 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,TTree* treev2){
889 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
890 //(or AliITStrackV2 tracks found with function FindTracks of this class)
895 if(fVert)delete fVert;
896 fVert = fVertexer->FindVertexForCurrentEvent(evnum);
899 gAlice->GetEvent(evnum);
901 Fatal("FindTracks","Vertex is missing\n");
905 Double_t primaryVertex[3];
906 fVert->GetXYZ(primaryVertex);
909 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
910 fTable->FillArray(fITSclusters);
911 fTable->FillArrayCoorAngles();
914 TBranch* bra = (TBranch*)treev2->GetBranch("tracks");
915 if(!bra) Warning("UseFoundTracksV2","No branch for track tree");
916 AliITStrackV2* ttrrt = new AliITStrackV2;
917 bra->SetAddress(&ttrrt);
919 for(Int_t nj=0;nj<treev2->GetEntries();nj++){
920 treev2->GetEvent(nj);
921 Int_t ncl = ttrrt->GetNumberOfClusters();
922 for(Int_t k=0;k<ncl;k++){
923 Int_t index = ttrrt->GetClusterIndex(k);
924 AliITSclusterV2* clui = (AliITSclusterV2*)GetCluster(index);
925 if(clui->IsUsed()==0) clui->Use();
932 //_______________________________________________________________________
933 void AliITStrackerSA::UseFoundTracksV2(AliESD *event){
934 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
938 Double_t primaryVertex[3];
939 event->GetVertex()->GetXYZ(primaryVertex);
942 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
943 fTable->FillArray(fITSclusters);
944 fTable->FillArrayCoorAngles();
947 Int_t ntracks = event->GetNumberOfTracks();
949 AliESDtrack *esd=event->GetTrack(ntracks);
950 if ((esd->GetStatus()&
951 AliESDtrack::kITSin|AliESDtrack::kTPCin)==0) continue;
953 Int_t ncl = esd->GetITSclusters(idx);
954 for(Int_t clu=0; clu<ncl; clu++){
955 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(idx[clu]);
956 if(cl->IsUsed()==0) cl->Use();
960 Info("UseFoundTracksV2","Clusters of tracks prolonged from TPC deleted");
965 //_______________________________________________________
966 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Int_t pflag,AliITSclusterTable* table){
967 //function used to to find the clusters associated to the track
969 AliITSlayer &lay = fgLayers[layer];
970 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
971 for(Int_t i=0;i<fGeom->GetNlayers();i++){
972 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
976 Float_t cx1,cx2,cy1,cy2;
977 FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
978 Int_t fun = FindIntersection(fCoef1,fCoef2,fCoef3,-(lay.GetR()*lay.GetR()),cx1,cy1,cx2,cy2);
983 Double_t fi1 =TMath::ATan2(cy1,cx1);
984 Double_t fi2 =TMath::ATan2(cy2,cx2);
985 fPhiEstimate = ChoosePoint(fi1,fi2,fPhic);
988 Double_t zed = TMath::Tan(fLambdac)*lay.GetR()+zvertex;
989 Double_t zed1 = TMath::Tan(fLambdac+lambdawindow)*lay.GetR()+zvertex;
990 Double_t zed2 = TMath::Tan(fLambdac-lambdawindow)*lay.GetR()+zvertex;
992 Double_t fi = fPhiEstimate;
993 Int_t nmod = lay.FindDetectorIndex(fi,zed);
998 nmod += firstmod[layer];
1000 Int_t nm[8]={0,0,0,0,0,0,0,0};
1001 nm[0] = lay.FindDetectorIndex(fi+phiwindow,zed);
1002 nm[1] = lay.FindDetectorIndex(fi-phiwindow,zed);
1003 nm[2] = lay.FindDetectorIndex(fi,zed1);
1004 nm[3] = lay.FindDetectorIndex(fi,zed2);
1005 nm[4] = lay.FindDetectorIndex(fi+phiwindow,zed1);
1006 nm[5] = lay.FindDetectorIndex(fi-phiwindow,zed1);
1007 nm[6] = lay.FindDetectorIndex(fi+phiwindow,zed2);
1008 nm[7] = lay.FindDetectorIndex(fi-phiwindow,zed2);
1012 TArrayI* array =(TArrayI*)table->GetListOfClusters(nmod);
1013 TArrayI* listc = new TArrayI(array->GetSize());
1014 for(Int_t i=0;i<array->GetSize();i++){
1015 Int_t in=(Int_t)array->At(i);
1016 listc->AddAt(in,nn);
1023 for(Int_t h=k+1;h<8;h++){
1036 for(Int_t ii=0;ii<8;ii++){
1037 if(nm[ii]!=value && nm[ii]!=nmod && nm[ii]>=0){
1038 TArrayI* ar =(TArrayI*)table->GetListOfClusters(nm[ii]+firstmod[layer]);
1039 listc->Set(listc->GetSize()+ar->GetSize());
1040 for(Int_t j=0;j<ar->GetSize();j++){
1041 Int_t in=(Int_t)ar->At(j);
1042 listc->AddAt(in,nn);
1050 for(Int_t i=0;i<listc->GetSize();i++){
1051 Int_t index = (Int_t)listc->At(i);
1052 AliITSclusterV2* cllay = lay.GetCluster(index);
1053 if(cllay==0) continue;
1054 if(cllay->IsUsed()==1) continue;
1055 Double_t phi = fTable->GetPhiCluster(layer,index);
1056 Double_t lambda= fTable->GetLambdaCluster(layer,index);
1058 if(TMath::Abs(fLambdac-lambda)<lambdawindow &&
1059 TMath::Abs(fPhiEstimate-phi)<phiwindow){
1062 if(trs->GetNumberOfClustersSA()==15){
1067 trs->AddClusterSA(layer,index);
1070 fPointc[0]=fTable->GetXCluster(layer,index);
1071 fPointc[1]=fTable->GetYCluster(layer,index);
1083 //________________________________________________________________
1084 void AliITStrackerSA::UpdatePoints(){
1085 //update of points for the estimation of the curvature
1087 fPoint1[0]=fPoint2[0];
1088 fPoint1[1]=fPoint2[1];
1089 fPoint2[0]=fPoint3[0];
1090 fPoint2[1]=fPoint3[1];
1091 fPoint3[0]=fPointc[0];
1092 fPoint3[1]=fPointc[1];
1099 //___________________________________________________________________
1100 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){
1102 //given (x,y) of three recpoints (in global coordinates)
1103 //returns the parameters a,b,c of circonference x*x + y*y +a*x + b*y +c
1105 Float_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1106 if(den==0) return 0;
1107 a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1108 b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1109 c = -x1*x1-y1*y1-a*x1-b*y1;
1112 //__________________________________________________________________________
1113 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){
1115 //Finds the intersection between the circonference of the track and the circonference centered in (0,0) represented by one layer
1116 //c2 is -rlayer*rlayer
1120 Float_t aA = (b1*b1)/(a1*a1)+1;
1121 Float_t bB = (-2*m*b1/(a1*a1));
1122 Float_t cC = c2+(m*m)/(a1*a1);
1123 if((bB*bB-4*aA*cC)<0) return 0;
1125 y1 = (-bB+TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA);
1126 y2 = (-bB-TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA);
1127 x1 = (c2-c1-b1*y1)/a1;
1128 x2 = (c2-c1-b1*y2)/a1;
1132 //____________________________________________________________________
1133 Double_t AliITStrackerSA::Curvature(Double_t x1,Double_t y1,Double_t
1134 x2,Double_t y2,Double_t x3,Double_t y3){
1136 //calculates the curvature of track
1137 Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1138 if(den==0) return 0;
1139 Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1140 Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1141 Double_t c = -x1*x1-y1*y1-a*x1-b*y1;
1144 if((a*a+b*b-4*c)<0) return 0;
1145 Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.;
1146 if(rad==0) return 0;
1148 if((x1>0 && y1>0 && x1<xc)) rad*=-1;
1149 if((x1<0 && y1>0 && x1<xc)) rad*=-1;
1150 // if((x1<0 && y1<0 && x1<xc)) rad*=-1;
1151 // if((x1>0 && y1<0 && x1<xc)) rad*=-1;
1156 //____________________________________________________________________
1157 Double_t AliITStrackerSA::ChoosePoint(Double_t p1, Double_t p2, Double_t pp){
1159 //Returns the point closest to pp
1161 Double_t diff1 = p1-pp;
1162 Double_t diff2 = p2-pp;
1164 if(TMath::Abs(diff1)<TMath::Abs(diff2)) fPhiEstimate=p1;
1165 else fPhiEstimate=p2;
1166 return fPhiEstimate;
1171 //_________________________________________________________________
1172 AliITStrackV2* AliITStrackerSA::FindTrackLowChiSquare(TObjArray* tracklist, Int_t dim) const {
1173 // returns track with lowes chi square
1175 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(0);
1178 if(dim==0) return 0;
1179 Double_t * chi2 = new Double_t[dim];
1180 Int_t * index = new Int_t[dim];
1181 for(Int_t i=0;i<dim;i++){
1182 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(i);
1183 chi2[i]=trk->GetChi2();
1187 Int_t w=0;Double_t value;
1190 for(Int_t j=w+1;j<dim;j++){
1191 if(chi2[w]<chi2[j]){
1203 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(index[dim-1]);
1210 //__________________________________________________________
1211 Int_t AliITStrackerSA::FindLabel(Int_t l1, Int_t l2, Int_t l3, Int_t l4, Int_t l5, Int_t l6){
1213 //function used to determine the track label
1215 Int_t lb[6] = {l1,l2,l3,l4,l5,l6};
1216 Int_t aa[6]={1,1,1,1,1,1};
1219 Int_t k=0;Int_t w=0;Int_t num=6;
1220 if(lb[5]==-1) num=5;
1224 for(Int_t i=k+1;i<num;i++){
1226 if(lb[k]==lb[i] && aa[k]!=0){
1237 for(Int_t j=0;j<6;j++){
1249 if(num==6) return lb[5];
1253 //_____________________________________________________________________________
1254 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,
1255 Int_t gl12, Int_t gl13, Int_t gl14,Int_t gl15, Int_t gl16, Int_t gl17, Int_t gl18, Int_t numberofpoints){
1258 //function used to assign label to the found track. If track is fake, the label is negative
1260 Int_t lb0[6] = {gl1,gl2,gl3,gl4,gl5,gl6};
1261 Int_t lb1[6] = {gl7,gl8,gl9,gl10,gl11,gl12};
1262 Int_t lb2[6] = {gl13,gl14,gl15,gl16,gl17,gl18};
1263 Int_t ll=FindLabel(lb0[0],lb0[1],lb0[2],lb0[3],lb0[4],lb0[5]);
1264 Int_t lflag=0;Int_t num=6;
1265 if(lb0[5]==-1 && lb1[5]==-1 && lb2[5]==-1) num=5;
1267 for(Int_t i=0;i<num;i++){
1268 if(lb0[i]==ll || lb1[i]==ll || lb2[i]==ll) lflag+=1;
1271 if(lflag>=numberofpoints) return ll;
1277 //_____________________________________________________________________________
1278 void AliITStrackerSA::SetWindowSizes(Int_t n, Double_t *phi, Double_t *lam){
1279 // Set sizes of the phi and lambda windows used for track finding
1281 if(phi){ // user defined values
1282 fPhiWin = new Double_t[fNloop];
1283 fLambdaWin = new Double_t[fNloop];
1284 for(Int_t k=0;k<fNloop;k++){
1286 fLambdaWin[k]=lam[k];
1289 else { // default values
1291 Double_t phid[32] = {0.002,0.003,0.004,0.0045,0.0047,
1292 0.005,0.0053,0.0055,
1293 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1294 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1295 0.01,0.0105,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1297 Double_t lambdad[32] = {0.002,0.003,0.004,0.0045,0.0047,
1298 0.005,0.0053,0.0055,
1299 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1300 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1301 0.01,0.015,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1308 fPhiWin = new Double_t[fNloop];
1309 fLambdaWin = new Double_t[fNloop];
1311 for(Int_t k=0;k<fNloop;k++){
1313 fLambdaWin[k]=lambdad[k];