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....
411 //______________________________________________________________________
412 Int_t AliITStrackerSA::FindTracks(AliESD* event){
414 // Track finder using the ESD object
418 Fatal("FindTracks","ITS cluster tree is not accessed - Abort!!!\n Please use method SetClusterTree to pass the pointer to the tree\n");
423 Double_t errorsprimvert[3];
424 Double_t primaryVertex[3];
425 event->GetVertex()->GetXYZ(primaryVertex);
426 event->GetVertex()->GetSigmaXYZ(errorsprimvert);
428 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
429 // Warning("FindTracks","Set errors on vertex positions x and y at 0.005");
430 errorsprimvert[0]=0.005;
431 errorsprimvert[1]=0.005;
434 //Fill array with cluster indices for each module
436 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
437 fTable->FillArray(fITSclusters);
438 fTable->FillArrayCoorAngles();
441 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
442 for(Int_t i=0;i<fGeom->GetNlayers();i++){
443 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
445 // firstmod [i] number of the first module in the ITS layer i.
448 AliITSlayer &layer=fgLayers[0]; // first layer
451 Int_t dim=layer.GetNumberOfClusters();
452 //loop on the different windows
453 for(Int_t nloop=0;nloop<fNloop;nloop++){
454 for(Int_t ncl=0;ncl<dim;ncl++){ //loop starting from layer 0
457 AliITSclusterV2* cl = (AliITSclusterV2*)layer.GetCluster(ncl);
459 if(cl->IsUsed()==1) continue;
460 fPhic = fTable->GetPhiCluster(0,ncl);
461 fLambdac = fTable->GetLambdaCluster(0,ncl);
462 fPhiEstimate = fPhic;
463 AliITStrackSA* trs = new AliITStrackSA();
464 fPoint1[0]=primaryVertex[0];
465 fPoint1[1]=primaryVertex[1];
466 fPoint2[0]=fTable->GetXCluster(0,ncl);
467 fPoint2[1]=fTable->GetYCluster(0,ncl);
468 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
469 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
470 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
472 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
475 fPoint3[0] = fPointc[0];
476 fPoint3[1] = fPointc[1];
478 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
479 if(nn[1]==0 && nn[2]==0) pflag=0;
480 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
481 if(nn[2]!=0 && nn[1]==0){
483 fPoint3[0]=fPointc[0];
484 fPoint3[1]=fPointc[1];
487 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
489 if(nn[3]!=0) UpdatePoints();
490 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
492 if(nn[4]!=0) UpdatePoints();
493 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
497 Int_t numberofpoints;
498 if(fSixPoints) numberofpoints=6; //check of the candidate track
499 else numberofpoints=5; //if track is good (with the required number
500 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
501 if(nn[nnp]!=0) layOK+=1;
503 if(layOK>=numberofpoints){
504 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
506 Int_t nct = trs->GetNumberOfClustersSA();
508 Int_t index = trs->GetClusterIndexSA(nct);
509 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
510 if(kl->IsUsed()==1) kl->Use();
515 AliESDtrack outtrack;
516 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
517 event->AddTrack(&outtrack);
519 Int_t nct = tr2->GetNumberOfClusters();
521 Int_t index = tr2->GetClusterIndex(nct);
522 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
523 if(kl->IsUsed()==0) kl->Use();
528 Int_t nct = trs->GetNumberOfClustersSA();
530 Int_t index = trs->GetClusterIndexSA(nct);
531 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
532 if(kl->IsUsed()==1) kl->Use();
538 }//end loop on clusters of layer1
542 //if 5/6 points are required, second loop starting
543 //from second layer, to find tracks with point of
547 // counter for clusters on each layer
548 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
549 for(Int_t nloop=0;nloop<fNloop;nloop++){
550 AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
551 Int_t ncl2=layer2.GetNumberOfClusters();
552 while(ncl2--){ //loop starting from layer 2
555 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
556 if(cl->IsUsed()==1) continue;
557 fPhic = fTable->GetPhiCluster(1,ncl2);
558 fLambdac = fTable->GetLambdaCluster(1,ncl2);
559 fPhiEstimate = fPhic;
560 AliITStrackSA* trs = new AliITStrackSA();
561 fPoint1[0]=primaryVertex[0];
562 fPoint1[1]=primaryVertex[1];
563 fPoint2[0]=fTable->GetXCluster(1,ncl2);
564 fPoint2[1]=fTable->GetYCluster(1,ncl2);
566 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
567 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
568 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
569 trs,primaryVertex[2],pflag,fTable);
577 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
578 if(nn[nnp]!=0) fl+=1;
581 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert);
584 Int_t nct = trs->GetNumberOfClustersSA();
586 Int_t index = trs->GetClusterIndexSA(nct);
587 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
588 if(kl->IsUsed()==1) kl->Use();
594 AliESDtrack outtrack;
595 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
596 event->AddTrack(&outtrack);
597 Int_t nct = tr2->GetNumberOfClusters();
599 Int_t index = tr2->GetClusterIndex(nct);
600 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
602 if(kl->IsUsed()==0) kl->Use();
606 Int_t nct = trs->GetNumberOfClustersSA();
608 Int_t index = trs->GetClusterIndexSA(nct);
609 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
611 if(kl->IsUsed()==1) kl->Use();
615 }//end loop on clusters of layer2
622 Info("FindTracks","Number of found tracks: %d",event->GetNumberOfTracks());
629 //________________________________________________________________________
631 AliITStrackV2* AliITStrackerSA::FitTrack(AliITStrackSA* tr,Double_t *primaryVertex,Double_t *errorsprimvert){
632 //fit of the found track
635 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
636 for(Int_t i=0;i<fGeom->GetNlayers();i++){
637 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
639 AliITStrackV2* otrack2;
640 Int_t nclusters = tr->GetNumberOfClustersSA();
641 TObjArray** listlayer = new TObjArray*[fGeom->GetNlayers()];
642 for(Int_t i=0;i<fGeom->GetNlayers();i++){
643 listlayer[i] = new TObjArray(0,0);
653 Int_t * nnn = new Int_t[fGeom->GetNlayers()];
654 for(Int_t i=0;i<fGeom->GetNlayers();i++)nnn[i]=0;
656 for(Int_t ncl=0;ncl<nclusters;ncl++){
657 Int_t index = tr->GetClusterIndexSA(ncl);
658 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(index);
660 if(cl->IsUsed()==1) cl->Use();
661 Int_t lay = (index & 0xf0000000) >> 28;
662 if(lay==0) { listlayer[0]->AddLast(cl); clind0[nnn[0]]=index;nnn[0]++;}
663 if(lay==1) { listlayer[1]->AddLast(cl); clind1[nnn[1]]=index;nnn[1]++;}
664 if(lay==2) { listlayer[2]->AddLast(cl); clind2[nnn[2]]=index;nnn[2]++;}
665 if(lay==3) { listlayer[3]->AddLast(cl); clind3[nnn[3]]=index;nnn[3]++;}
666 if(lay==4) { listlayer[4]->AddLast(cl); clind4[nnn[4]]=index;nnn[4]++;}
667 if(lay==5) { listlayer[5]->AddLast(cl); clind5[nnn[5]]=index;nnn[5]++;}
672 Int_t * end = new Int_t[fGeom->GetNlayers()];
673 for(Int_t i=0;i<fGeom->GetNlayers();i++){
674 if(listlayer[i]->GetEntries()==0) end[i]=1;
675 else end[i]=listlayer[i]->GetEntries();
678 for(Int_t l1=0;l1<end[0];l1++){//loop on layer 1
679 AliITSclusterV2* cl0 = (AliITSclusterV2*)listlayer[0]->At(l1);
680 TVector3** recp = new TVector3*[3];
681 TVector3** errs = new TVector3*[3];
682 recp[0] = new TVector3(primaryVertex[0],primaryVertex[1],primaryVertex[2]);
683 errs[0] = new TVector3(errorsprimvert[0],errorsprimvert[1],errorsprimvert[2]);
684 Double_t x1,y1,z1,sx1,sy1,sz1;
685 Double_t x2,y2,z2,sx2,sy2,sz2;
686 AliITSclusterV2* p1=0;
687 AliITSclusterV2* p2=0;
688 Int_t index1=clind0[l1];
690 for(Int_t l2=0;l2<end[1];l2++){//loop on layer 2
691 AliITSclusterV2* cl1 = (AliITSclusterV2*)listlayer[1]->At(l2);
693 for(Int_t l3=0;l3<end[2];l3++){ //loop on layer 3
694 AliITSclusterV2* cl2 = (AliITSclusterV2*)listlayer[2]->At(l3);
696 if(cl0==0 && cl1!=0) {
697 p2 = cl2;index1=clind2[l3];
701 if(cl0!=0 && cl1==0){
703 p2=cl2;index2=clind2[l3];
705 if(cl0!=0 && cl1!=0){
709 Int_t lay1=(index1 & 0xf0000000) >> 28;
710 Int_t cln1=(index1 & 0x0fffffff) >> 00;
711 Int_t lay2=(index2 & 0xf0000000) >> 28;
712 Int_t cln2=(index2 & 0x0fffffff) >> 00;
713 x1 = fTable->GetXCluster(lay1,cln1);
714 x2 = fTable->GetXCluster(lay2,cln2);
715 y1 = fTable->GetYCluster(lay1,cln1);
716 y2 = fTable->GetYCluster(lay2,cln2);
717 z1 = fTable->GetZCluster(lay1,cln1);
718 z2 = fTable->GetZCluster(lay2,cln2);
719 sx1 = fTable->GetXClusterError(lay1,cln1);
720 sx2 = fTable->GetXClusterError(lay2,cln2);
721 sy1 = fTable->GetYClusterError(lay1,cln1);
722 sy2 = fTable->GetYClusterError(lay2,cln2);
723 sz1 = fTable->GetZClusterError(lay1,cln1);
724 sz2 = fTable->GetZClusterError(lay2,cln2);
725 Double_t phi1 = fTable->GetPhiCluster(lay1,cln1);
726 Int_t module1 = p1->GetDetectorIndex()+firstmod[0];
727 recp[1] = new TVector3(x1,y1,z1);
728 errs[1] = new TVector3(sx1,sy1,sz1);
729 recp[2] = new TVector3(x2,y2,z2);
730 errs[2] = new TVector3(sx2,sy2,sz2);
732 //fit on the Riemann sphere
733 Float_t seed1,seed2,seed3;
734 AliITSRiemannFit fit;
735 Int_t rf = fit.FitHelix(3,recp,errs,seed1,seed2,seed3); //this gives phi,tgl,curvature to start Kalman Filter
741 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
742 phi=seed1+1.5*TMath::Pi();
745 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
746 phi=seed1+0.5*TMath::Pi();
749 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
750 phi=seed1-0.5*TMath::Pi();
755 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
756 phi=seed1+0.5*TMath::Pi();
759 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
760 phi=seed1-0.5*TMath::Pi();
763 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
764 phi=seed1-1.5*TMath::Pi();
769 Int_t layer,ladder,detector;
770 fGeom->GetModuleId(module1,layer,ladder,detector);
771 Float_t yclu1 = p1->GetY();
772 Float_t zclu1 = p1->GetZ();
773 Double_t cv=Curvature(primaryVertex[0],primaryVertex[1],x1,y1,x2,y2);
775 for(Int_t l4=0;l4<end[3];l4++){ //loop on layer 4
776 AliITSclusterV2* cl3 = (AliITSclusterV2*)listlayer[3]->At(l4);
777 for(Int_t l5=0;l5<end[4];l5++){ //loop on layer 5
778 AliITSclusterV2* cl4 = (AliITSclusterV2*)listlayer[4]->At(l5);
779 for(Int_t l6=0;l6<end[5];l6++){ //loop on layer 6
780 AliITSclusterV2* cl5 = (AliITSclusterV2*)listlayer[5]->At(l6);
781 AliITStrackSA* trac = new AliITStrackSA(layer,ladder,detector,yclu1,zclu1,phi,tgl,cv,1);
783 if(cl5!=0) trac->AddClusterV2(5,(clind5[l6] & 0x0fffffff)>>0);
784 if(cl4!=0) trac->AddClusterV2(4,(clind4[l5] & 0x0fffffff)>>0);
785 if(cl3!=0) trac->AddClusterV2(3,(clind3[l4] & 0x0fffffff)>>0);
786 if(cl2!=0) trac->AddClusterV2(2,(clind2[l3] & 0x0fffffff)>>0);
787 if(cl1!=0) trac->AddClusterV2(1,(clind1[l2] & 0x0fffffff)>>0);
788 if(cl0!=0) trac->AddClusterV2(0,(clind0[l1] & 0x0fffffff)>>0);
790 //fit with Kalman filter using AliITStrackerV2::RefitAt()
792 AliITStrackV2* ot = new AliITStrackV2(*trac);
794 ot->ResetCovariance();
797 if(RefitAt(49.,ot,trac)){ //fit from layer 1 to layer 6
799 otrack2 = new AliITStrackV2(*ot);
800 otrack2->ResetCovariance();
801 otrack2->ResetClusters();
802 //fit from layer 6 to layer 1
803 if(RefitAt(3.7,otrack2,ot)) fListOfTracks->AddLast(otrack2);
813 for(Int_t i=1;i<3;i++){
827 Int_t dim=fListOfTracks->GetEntries();
829 for(Int_t i=0;i<fGeom->GetNlayers();i++){
836 AliITStrackV2* otrack =(AliITStrackV2*)FindTrackLowChiSquare(fListOfTracks,dim);
838 if(otrack==0) return 0;
839 Int_t * indexc = new Int_t[fGeom->GetNlayers()];
840 for(Int_t i=0;i<fGeom->GetNlayers();i++) indexc[i]=0;
841 for(Int_t nind=0;nind<otrack->GetNumberOfClusters();nind++){
842 indexc[nind] = otrack->GetClusterIndex(nind);
844 AliITSclusterV2* cl0 = (AliITSclusterV2*)GetCluster(indexc[0]);
845 AliITSclusterV2* cl1 = (AliITSclusterV2*)GetCluster(indexc[1]);
846 AliITSclusterV2* cl2 = (AliITSclusterV2*)GetCluster(indexc[2]);
847 AliITSclusterV2* cl3 = (AliITSclusterV2*)GetCluster(indexc[3]);
848 AliITSclusterV2* cl4 = (AliITSclusterV2*)GetCluster(indexc[4]);
849 Int_t labl[3]={-1,-1,-1};
850 if(otrack->GetNumberOfClusters()==fGeom->GetNlayers()){
851 AliITSclusterV2* cl5 = (AliITSclusterV2*)GetCluster(indexc[5]);
852 labl[0]=cl5->GetLabel(0);
853 labl[1]=cl5->GetLabel(1);
854 labl[2]=cl5->GetLabel(2);
857 if(otrack->GetNumberOfClusters()==(fGeom->GetNlayers()-1)){
862 Int_t numberofpoints;
863 if(fSixPoints) numberofpoints=6;
864 else numberofpoints=5;
865 Int_t label = Label(cl0->GetLabel(0),cl1->GetLabel(0),
866 cl2->GetLabel(0),cl3->GetLabel(0),
867 cl4->GetLabel(0),labl[0],
868 cl0->GetLabel(1),cl1->GetLabel(1),
869 cl2->GetLabel(1),cl3->GetLabel(1),
870 cl4->GetLabel(1),labl[1],
871 cl0->GetLabel(2),cl1->GetLabel(2),
872 cl2->GetLabel(2),cl3->GetLabel(2),
873 cl4->GetLabel(2),labl[2],numberofpoints);
875 otrack->SetLabel(label);
876 for(Int_t i=0;i<fGeom->GetNlayers();i++){
885 //_______________________________________________________________________
886 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,TTree* treev2){
887 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
888 //(or AliITStrackV2 tracks found with function FindTracks of this class)
893 if(fVert)delete fVert;
894 fVert = fVertexer->FindVertexForCurrentEvent(evnum);
897 gAlice->GetEvent(evnum);
899 Fatal("FindTracks","Vertex is missing\n");
903 Double_t primaryVertex[3];
904 fVert->GetXYZ(primaryVertex);
907 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
908 fTable->FillArray(fITSclusters);
909 fTable->FillArrayCoorAngles();
912 TBranch* bra = (TBranch*)treev2->GetBranch("tracks");
913 if(!bra) Warning("UseFoundTracksV2","No branch for track tree");
914 AliITStrackV2* ttrrt = new AliITStrackV2;
915 bra->SetAddress(&ttrrt);
917 for(Int_t nj=0;nj<treev2->GetEntries();nj++){
918 treev2->GetEvent(nj);
919 Int_t ncl = ttrrt->GetNumberOfClusters();
920 for(Int_t k=0;k<ncl;k++){
921 Int_t index = ttrrt->GetClusterIndex(k);
922 AliITSclusterV2* clui = (AliITSclusterV2*)GetCluster(index);
923 if(clui->IsUsed()==0) clui->Use();
930 //_______________________________________________________________________
931 void AliITStrackerSA::UseFoundTracksV2(AliESD *event){
932 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
936 Double_t primaryVertex[3];
937 event->GetVertex()->GetXYZ(primaryVertex);
940 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
941 fTable->FillArray(fITSclusters);
942 fTable->FillArrayCoorAngles();
945 Int_t ntracks = event->GetNumberOfTracks();
947 AliESDtrack *esd=event->GetTrack(ntracks);
948 if ((esd->GetStatus()&
949 AliESDtrack::kITSin|AliESDtrack::kTPCin)==0) continue;
951 Int_t ncl = esd->GetITSclusters(idx);
952 for(Int_t clu=0; clu<ncl; clu++){
953 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(idx[clu]);
954 if(cl->IsUsed()==0) cl->Use();
958 Info("UseFoundTracksV2","Clusters of tracks prolonged from TPC deleted");
963 //_______________________________________________________
964 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Int_t pflag,AliITSclusterTable* table){
965 //function used to to find the clusters associated to the track
967 AliITSlayer &lay = fgLayers[layer];
968 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
969 for(Int_t i=0;i<fGeom->GetNlayers();i++){
970 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
974 Float_t cx1,cx2,cy1,cy2;
975 FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
976 Int_t fun = FindIntersection(fCoef1,fCoef2,fCoef3,-(lay.GetR()*lay.GetR()),cx1,cy1,cx2,cy2);
981 Double_t fi1 =TMath::ATan2(cy1,cx1);
982 Double_t fi2 =TMath::ATan2(cy2,cx2);
983 fPhiEstimate = ChoosePoint(fi1,fi2,fPhic);
986 Double_t zed = TMath::Tan(fLambdac)*lay.GetR()+zvertex;
987 Double_t zed1 = TMath::Tan(fLambdac+lambdawindow)*lay.GetR()+zvertex;
988 Double_t zed2 = TMath::Tan(fLambdac-lambdawindow)*lay.GetR()+zvertex;
990 Double_t fi = fPhiEstimate;
991 Int_t nmod = lay.FindDetectorIndex(fi,zed);
996 nmod += firstmod[layer];
998 Int_t nm[8]={0,0,0,0,0,0,0,0};
999 nm[0] = lay.FindDetectorIndex(fi+phiwindow,zed);
1000 nm[1] = lay.FindDetectorIndex(fi-phiwindow,zed);
1001 nm[2] = lay.FindDetectorIndex(fi,zed1);
1002 nm[3] = lay.FindDetectorIndex(fi,zed2);
1003 nm[4] = lay.FindDetectorIndex(fi+phiwindow,zed1);
1004 nm[5] = lay.FindDetectorIndex(fi-phiwindow,zed1);
1005 nm[6] = lay.FindDetectorIndex(fi+phiwindow,zed2);
1006 nm[7] = lay.FindDetectorIndex(fi-phiwindow,zed2);
1010 TArrayI* array =(TArrayI*)table->GetListOfClusters(nmod);
1011 TArrayI* listc = new TArrayI(array->GetSize());
1012 for(Int_t i=0;i<array->GetSize();i++){
1013 Int_t in=(Int_t)array->At(i);
1014 listc->AddAt(in,nn);
1021 for(Int_t h=k+1;h<8;h++){
1034 for(Int_t ii=0;ii<8;ii++){
1035 if(nm[ii]!=value && nm[ii]!=nmod && nm[ii]>=0){
1036 TArrayI* ar =(TArrayI*)table->GetListOfClusters(nm[ii]+firstmod[layer]);
1037 listc->Set(listc->GetSize()+ar->GetSize());
1038 for(Int_t j=0;j<ar->GetSize();j++){
1039 Int_t in=(Int_t)ar->At(j);
1040 listc->AddAt(in,nn);
1048 for(Int_t i=0;i<listc->GetSize();i++){
1049 Int_t index = (Int_t)listc->At(i);
1050 AliITSclusterV2* cllay = lay.GetCluster(index);
1051 if(cllay==0) continue;
1052 if(cllay->IsUsed()==1) continue;
1053 Double_t phi = fTable->GetPhiCluster(layer,index);
1054 Double_t lambda= fTable->GetLambdaCluster(layer,index);
1056 if(TMath::Abs(fLambdac-lambda)<lambdawindow &&
1057 TMath::Abs(fPhiEstimate-phi)<phiwindow){
1060 if(trs->GetNumberOfClustersSA()==15){
1065 trs->AddClusterSA(layer,index);
1068 fPointc[0]=fTable->GetXCluster(layer,index);
1069 fPointc[1]=fTable->GetYCluster(layer,index);
1081 //________________________________________________________________
1082 void AliITStrackerSA::UpdatePoints(){
1083 //update of points for the estimation of the curvature
1085 fPoint1[0]=fPoint2[0];
1086 fPoint1[1]=fPoint2[1];
1087 fPoint2[0]=fPoint3[0];
1088 fPoint2[1]=fPoint3[1];
1089 fPoint3[0]=fPointc[0];
1090 fPoint3[1]=fPointc[1];
1097 //___________________________________________________________________
1098 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){
1100 //given (x,y) of three recpoints (in global coordinates)
1101 //returns the parameters a,b,c of circonference x*x + y*y +a*x + b*y +c
1103 Float_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1104 if(den==0) return 0;
1105 a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1106 b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1107 c = -x1*x1-y1*y1-a*x1-b*y1;
1110 //__________________________________________________________________________
1111 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){
1113 //Finds the intersection between the circonference of the track and the circonference centered in (0,0) represented by one layer
1114 //c2 is -rlayer*rlayer
1118 Float_t aA = (b1*b1)/(a1*a1)+1;
1119 Float_t bB = (-2*m*b1/(a1*a1));
1120 Float_t cC = c2+(m*m)/(a1*a1);
1121 if((bB*bB-4*aA*cC)<0) return 0;
1123 y1 = (-bB+TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA);
1124 y2 = (-bB-TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA);
1125 x1 = (c2-c1-b1*y1)/a1;
1126 x2 = (c2-c1-b1*y2)/a1;
1130 //____________________________________________________________________
1131 Double_t AliITStrackerSA::Curvature(Double_t x1,Double_t y1,Double_t
1132 x2,Double_t y2,Double_t x3,Double_t y3){
1134 //calculates the curvature of track
1135 Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1136 if(den==0) return 0;
1137 Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1138 Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1139 Double_t c = -x1*x1-y1*y1-a*x1-b*y1;
1142 if((a*a+b*b-4*c)<0) return 0;
1143 Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.;
1144 if(rad==0) return 0;
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;
1149 // if((x1>0 && y1<0 && x1<xc)) rad*=-1;
1154 //____________________________________________________________________
1155 Double_t AliITStrackerSA::ChoosePoint(Double_t p1, Double_t p2, Double_t pp){
1157 //Returns the point closest to pp
1159 Double_t diff1 = p1-pp;
1160 Double_t diff2 = p2-pp;
1162 if(TMath::Abs(diff1)<TMath::Abs(diff2)) fPhiEstimate=p1;
1163 else fPhiEstimate=p2;
1164 return fPhiEstimate;
1169 //_________________________________________________________________
1170 AliITStrackV2* AliITStrackerSA::FindTrackLowChiSquare(TObjArray* tracklist, Int_t dim) const {
1171 // returns track with lowes chi square
1173 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(0);
1176 if(dim==0) return 0;
1177 Double_t * chi2 = new Double_t[dim];
1178 Int_t * index = new Int_t[dim];
1179 for(Int_t i=0;i<dim;i++){
1180 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(i);
1181 chi2[i]=trk->GetChi2();
1185 Int_t w=0;Double_t value;
1188 for(Int_t j=w+1;j<dim;j++){
1189 if(chi2[w]<chi2[j]){
1201 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(index[dim-1]);
1208 //__________________________________________________________
1209 Int_t AliITStrackerSA::FindLabel(Int_t l1, Int_t l2, Int_t l3, Int_t l4, Int_t l5, Int_t l6){
1211 //function used to determine the track label
1213 Int_t lb[6] = {l1,l2,l3,l4,l5,l6};
1214 Int_t aa[6]={1,1,1,1,1,1};
1217 Int_t k=0;Int_t w=0;Int_t num=6;
1218 if(lb[5]==-1) num=5;
1222 for(Int_t i=k+1;i<num;i++){
1224 if(lb[k]==lb[i] && aa[k]!=0){
1235 for(Int_t j=0;j<6;j++){
1247 if(num==6) return lb[5];
1251 //_____________________________________________________________________________
1252 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,
1253 Int_t gl12, Int_t gl13, Int_t gl14,Int_t gl15, Int_t gl16, Int_t gl17, Int_t gl18, Int_t numberofpoints){
1256 //function used to assign label to the found track. If track is fake, the label is negative
1258 Int_t lb0[6] = {gl1,gl2,gl3,gl4,gl5,gl6};
1259 Int_t lb1[6] = {gl7,gl8,gl9,gl10,gl11,gl12};
1260 Int_t lb2[6] = {gl13,gl14,gl15,gl16,gl17,gl18};
1261 Int_t ll=FindLabel(lb0[0],lb0[1],lb0[2],lb0[3],lb0[4],lb0[5]);
1262 Int_t lflag=0;Int_t num=6;
1263 if(lb0[5]==-1 && lb1[5]==-1 && lb2[5]==-1) num=5;
1265 for(Int_t i=0;i<num;i++){
1266 if(lb0[i]==ll || lb1[i]==ll || lb2[i]==ll) lflag+=1;
1269 if(lflag>=numberofpoints) return ll;
1275 //_____________________________________________________________________________
1276 void AliITStrackerSA::SetWindowSizes(Int_t n, Double_t *phi, Double_t *lam){
1277 // Set sizes of the phi and lambda windows used for track finding
1279 if(phi){ // user defined values
1280 fPhiWin = new Double_t[fNloop];
1281 fLambdaWin = new Double_t[fNloop];
1282 for(Int_t k=0;k<fNloop;k++){
1284 fLambdaWin[k]=lam[k];
1287 else { // default values
1289 Double_t phid[32] = {0.002,0.003,0.004,0.0045,0.0047,
1290 0.005,0.0053,0.0055,
1291 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1292 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1293 0.01,0.0105,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1295 Double_t lambdad[32] = {0.002,0.003,0.004,0.0045,0.0047,
1296 0.005,0.0053,0.0055,
1297 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1298 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1299 0.01,0.015,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1306 fPhiWin = new Double_t[fNloop];
1307 fLambdaWin = new Double_t[fNloop];
1309 for(Int_t k=0;k<fNloop;k++){
1311 fLambdaWin[k]=lambdad[k];