1 /**************************************************************************
2 * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 ////////////////////////////////////////////////////
17 // Stand alone tracker class //
18 // Origin: Elisabetta Crescio //
19 // e-mail: crescio@to.infn.it //
20 // tracks are saved as AliITStrackV2 objects //
21 ////////////////////////////////////////////////////
27 #include <TObjArray.h>
30 #include "AliITSclusterTable.h"
31 #include "AliITSclusterV2.h"
32 #include "AliITSgeom.h"
33 #include "AliITSRiemannFit.h"
34 #include "AliITStrackerSA.h"
35 #include "AliITStrackSA.h"
36 #include "AliITSVertexer.h"
37 #include "AliITSVertex.h"
39 #include "AliESDtrack.h"
41 ClassImp(AliITStrackerSA)
43 //____________________________________________________________________________
44 AliITStrackerSA::AliITStrackerSA():AliITStrackerV2(){
45 // Default constructor
48 //____________________________________________________________________________
49 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertex *vert):AliITStrackerV2(geom)
51 // Standard constructor (Vertex is known and passed to this obj.)
57 //______________________________________________________________________
58 AliITStrackerSA::AliITStrackerSA(const AliITStrackerSA &trkr) :
59 AliITStrackerV2(trkr) {
61 // Copies are not allowed. The method is protected to avoid misuse.
62 Error("AliITStrackerSA","Copy constructor not allowed\n");
65 //______________________________________________________________________
66 AliITStrackerSA& AliITStrackerSA::operator=(const
67 AliITStrackerSA& /* trkr */){
68 // Assignment operator
69 // Assignment is not allowed. The method is protected to avoid misuse.
70 Error("= operator","Assignment operator not allowed\n");
74 //____________________________________________________________________________
75 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertexer *vertexer):AliITStrackerV2(geom)
77 // Standard constructor (Vertex is unknown - vertexer is passed to this obj)
83 //____________________________________________________________________________
84 AliITStrackerSA::AliITStrackerSA(AliITStrackerSA& tracker):AliITStrackerV2(){
86 fPhiEstimate = tracker.fPhiEstimate;
87 for(Int_t i=0;i<2;i++){
88 fPoint1[i]=tracker.fPoint1[i];
89 fPoint2[i]=tracker.fPoint2[i];
90 fPoint3[i]=tracker.fPoint3[i];
91 fPointc[i]=tracker.fPointc[i];
93 fLambdac = tracker.fLambdac;
94 fPhic = tracker.fPhic;
95 fCoef1 = tracker.fCoef1;
96 fCoef2 = tracker.fCoef2;
97 fCoef3 = tracker.fCoef3;
98 fNloop = tracker.fNloop;
99 fPhiWin = tracker.fPhiWin;
100 fLambdaWin = tracker.fLambdaWin;
101 if(tracker.fVertexer && tracker.fVert){
102 fVert = new AliITSVertex(*tracker.fVert);
105 fVert = tracker.fVert;
107 fVertexer = tracker.fVertexer;
108 fGeom = tracker.fGeom;
109 fFlagLoad = tracker.fFlagLoad;
110 fTable = tracker.fTable;
113 //____________________________________________________________________________
114 AliITStrackerSA::~AliITStrackerSA(){
116 // if fVertexer is not null, the AliITSVertex obj. is owned by this class
117 // and is deleted here
119 if(fVert)delete fVert;
125 if(fPhiWin)delete []fPhiWin;
126 if(fLambdaWin)delete []fLambdaWin;
130 //____________________________________________________________________________
131 void AliITStrackerSA::Init(){
132 // Reset all data members
134 for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
149 //_______________________________________________________________________
150 void AliITStrackerSA::ResetForFinding(){
151 // Reset data members used in all loops during track finding
153 for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
162 //____________________________________________________________________________
163 void AliITStrackerSA::FindTracks(TTree *clusterTree, TTree *out,Int_t evnumber,char *opt){
165 /**************************************************************************
166 * This function finds primary tracks.
167 * Options: to use only ITS execute FindTracks *
168 * to define good tracks with 6 points out of 6 in the ITS write *
169 * *opt="6/6". The default is *opt="6/6" *
172 * Example: to execute function with only the ITS (no combined tracking *
173 * with TPC+ITS) and requiring 5/6 points to define a good track *
174 * use: FindTracks(treein,treeout,evnumber,"5/6") *
175 * to execute combined tracking, before using FindTracks use *
177 *************************************************************************/
181 Bool_t sixpoints= choice.Contains("6/6");
185 if(fVert)delete fVert;
186 fVert = fVertexer->FindVertexForCurrentEvent(evnumber);
189 gAlice->GetEvent(evnumber);
191 Fatal("FindTracks","Vertex is missing\n");
195 Double_t primaryVertex[3];
196 Double_t errorsprimvert[3];
197 fVert->GetXYZ(primaryVertex);
198 fVert->GetSigmaXYZ(errorsprimvert);
199 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
200 Warning("FindTracks","Set errors on vertex positions x and y at 0.0001");
201 errorsprimvert[0]=0.0001;
202 errorsprimvert[1]=0.0001;
204 fVert->PrintStatus();
207 //Fill array with cluster indices for each module
209 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
210 fTable->FillArray(clusterTree,evnumber);
211 fTable->FillArrayCoorAngles(clusterTree,evnumber);
214 SetEventNumber(evnumber);
215 if(GetFlagLoadedClusters()==0) LoadClusters(clusterTree);
217 //Fill tree for found tracks
218 AliITStrackV2* outrack=0;
219 TBranch* branch=out->Branch("tracks","AliITStrackV2",&outrack,32000,0);
220 if (!branch) out->Branch("tracks","AliITStrackV2",&outrack,32000,3);
221 else branch->SetAddress(&outrack);
224 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
225 for(Int_t i=0;i<fGeom->GetNlayers();i++){
226 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
228 // firstmod [i] number of the first module in the ITS layer i.
230 AliITSlayer &layer=fgLayers[0]; // first layer
232 Int_t dim=layer.GetNumberOfClusters();
233 //loop on the different windows
234 for(Int_t nloop=0;nloop<fNloop;nloop++){
235 for(Int_t ncl=0;ncl<dim;ncl++){
236 //loop starting from layer 0
239 AliITSclusterV2* cl = layer.GetCluster(ncl);
240 if(cl->IsUsed()==1) continue;
241 fPhic = fTable->GetPhiCluster(0,ncl);
242 fLambdac = fTable->GetLambdaCluster(0,ncl);
243 fPhiEstimate = fPhic;
244 AliITStrackSA* trs = new AliITStrackSA();
245 fPoint1[0]=primaryVertex[0];
246 fPoint1[1]=primaryVertex[1];
247 fPoint2[0]=fTable->GetXCluster(0,ncl);;
248 fPoint2[1]=fTable->GetYCluster(0,ncl);;
250 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
251 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
252 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
253 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
257 fPoint3[0] = fPointc[0];
258 fPoint3[1] = fPointc[1];
260 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
261 if(nn[1]==0 && nn[2]==0) pflag=0;
262 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
263 if(nn[2]!=0 && nn[1]==0){
265 fPoint3[0]=fPointc[0];
266 fPoint3[1]=fPointc[1];
269 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
271 if(nn[3]!=0) UpdatePoints();
272 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
274 if(nn[4]!=0) UpdatePoints();
275 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
279 Int_t numberofpoints;
280 if(sixpoints) numberofpoints=6; //check of the candidate track
281 else numberofpoints=5; //if track is good (with the required number
282 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
283 if(nn[nnp]!=0) layOK+=1;
285 if(layOK>=numberofpoints){
286 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert,opt);
288 Int_t nct = trs->GetNumberOfClustersSA();
290 Int_t index = trs->GetClusterIndexSA(nct);
291 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
292 if(kl->IsUsed()==1) kl->Use();
299 Int_t nct = tr2->GetNumberOfClusters();
302 Int_t index = tr2->GetClusterIndex(nct);
303 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
304 if(kl->IsUsed()==0) kl->Use();
309 Int_t nct = trs->GetNumberOfClustersSA();
311 Int_t index = trs->GetClusterIndexSA(nct);
312 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
313 if(kl->IsUsed()==1) kl->Use();
318 }//end loop on clusters of layer1
322 //if 5/6 points are required, second loop starting
323 //from second layer, to find tracks with point of
327 // counter for clusters on each layer
328 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
329 for(Int_t nloop=0;nloop<fNloop;nloop++){
330 AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
331 Int_t ncl2=layer2.GetNumberOfClusters();
332 while(ncl2--){ //loop starting from layer 2
335 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
336 if(cl->IsUsed()==1) continue;
337 fPhic = fTable->GetPhiCluster(1,ncl2);
338 fLambdac = fTable->GetLambdaCluster(1,ncl2);
339 fPhiEstimate = fPhic;
340 AliITStrackSA* trs = new AliITStrackSA();
342 fPoint1[0]=primaryVertex[0];
343 fPoint1[1]=primaryVertex[1];
344 fPoint2[0]=fTable->GetXCluster(1,ncl2);;
345 fPoint2[1]=fTable->GetYCluster(1,ncl2);;
347 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
348 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
349 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
350 trs,primaryVertex[2],pflag,fTable);
358 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
359 if(nn[nnp]!=0) fl+=1;
362 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert,opt);
364 Int_t nct = trs->GetNumberOfClustersSA();
366 Int_t index = trs->GetClusterIndexSA(nct);
367 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
368 if(kl->IsUsed()==1) kl->Use();
374 Int_t nct = tr2->GetNumberOfClusters();
376 Int_t index = tr2->GetClusterIndex(nct);
377 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
379 if(kl->IsUsed()==0) kl->Use();
383 Int_t nct = trs->GetNumberOfClustersSA();
385 Int_t index = trs->GetClusterIndexSA(nct);
386 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
388 if(kl->IsUsed()==1) kl->Use();
392 }//end loop on clusters of layer2
402 //______________________________________________________________________
403 void AliITStrackerSA::FindTracks(TTree *clusterTree, AliESD* event, Int_t evnumber,char *opt){
404 /**************************************************************************
405 * Options: to use only ITS execute FindTracks *
406 * to define good tracks with 6 points out of 6 in the ITS write *
407 * *opt="6/6". The default is *opt="6/6" *
410 * Example: to execute function with only the ITS (no combined tracking *
411 * with TPC+ITS) and requiring 5/6 points to define a good track *
412 * use: FindTracks(treein,treeout,evnumber,"5/6") *
413 * to execute combined tracking, before using FindTracks use *
415 *************************************************************************/
419 Bool_t sixpoints= choice.Contains("6/6");
422 if(fVert)delete fVert;
423 fVert = fVertexer->FindVertexForCurrentEvent(evnumber);
426 gAlice->GetEvent(evnumber);
428 Fatal("FindTracks","Vertex is missing\n");
432 Double_t primaryVertex[3];
433 Double_t errorsprimvert[3];
434 fVert->GetXYZ(primaryVertex);
435 fVert->GetSigmaXYZ(errorsprimvert);
436 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
437 Warning("FindTracks","Set errors on vertex positions x and y at 0.0001");
438 errorsprimvert[0]=0.0001;
439 errorsprimvert[1]=0.0001;
441 fVert->PrintStatus();
443 //Fill array with cluster indices for each module
445 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
446 fTable->FillArray(clusterTree,evnumber);
447 fTable->FillArrayCoorAngles(clusterTree,evnumber);
450 SetEventNumber(evnumber);
451 if(GetFlagLoadedClusters()==0) LoadClusters(clusterTree);
453 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
454 for(Int_t i=0;i<fGeom->GetNlayers();i++){
455 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
457 // firstmod [i] number of the first module in the ITS layer i.
460 AliITSlayer &layer=fgLayers[0]; // first layer
463 Int_t dim=layer.GetNumberOfClusters();
464 //loop on the different windows
465 for(Int_t nloop=0;nloop<fNloop;nloop++){
466 for(Int_t ncl=0;ncl<dim;ncl++){ //loop starting from layer 0
469 AliITSclusterV2* cl = (AliITSclusterV2*)layer.GetCluster(ncl);
471 if(cl->IsUsed()==1) continue;
472 fPhic = fTable->GetPhiCluster(0,ncl);
473 fLambdac = fTable->GetLambdaCluster(0,ncl);
474 fPhiEstimate = fPhic;
475 AliITStrackSA* trs = new AliITStrackSA();
476 fPoint1[0]=primaryVertex[0];
477 fPoint1[1]=primaryVertex[1];
478 fPoint2[0]=fTable->GetXCluster(0,ncl);
479 fPoint2[1]=fTable->GetYCluster(0,ncl);
480 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
481 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
482 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
484 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
487 fPoint3[0] = fPointc[0];
488 fPoint3[1] = fPointc[1];
490 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
491 if(nn[1]==0 && nn[2]==0) pflag=0;
492 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
493 if(nn[2]!=0 && nn[1]==0){
495 fPoint3[0]=fPointc[0];
496 fPoint3[1]=fPointc[1];
499 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
501 if(nn[3]!=0) UpdatePoints();
502 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
504 if(nn[4]!=0) UpdatePoints();
505 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],pflag,fTable);
509 Int_t numberofpoints;
510 if(sixpoints) numberofpoints=6; //check of the candidate track
511 else numberofpoints=5; //if track is good (with the required number
512 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
513 if(nn[nnp]!=0) layOK+=1;
515 if(layOK>=numberofpoints){
516 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert,opt);
518 Int_t nct = trs->GetNumberOfClustersSA();
520 Int_t index = trs->GetClusterIndexSA(nct);
521 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
522 if(kl->IsUsed()==1) kl->Use();
527 AliESDtrack outtrack;
528 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
529 event->AddTrack(&outtrack);
532 Int_t nct = tr2->GetNumberOfClusters();
534 Int_t index = tr2->GetClusterIndex(nct);
535 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
536 if(kl->IsUsed()==0) kl->Use();
540 Int_t nct = trs->GetNumberOfClustersSA();
542 Int_t index = trs->GetClusterIndexSA(nct);
543 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
544 if(kl->IsUsed()==1) kl->Use();
549 }//end loop on clusters of layer1
553 //if 5/6 points are required, second loop starting
554 //from second layer, to find tracks with point of
558 // counter for clusters on each layer
559 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
560 for(Int_t nloop=0;nloop<fNloop;nloop++){
561 AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
562 Int_t ncl2=layer2.GetNumberOfClusters();
563 while(ncl2--){ //loop starting from layer 2
566 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
567 if(cl->IsUsed()==1) continue;
568 fPhic = fTable->GetPhiCluster(1,ncl2);
569 fLambdac = fTable->GetLambdaCluster(1,ncl2);
570 fPhiEstimate = fPhic;
571 AliITStrackSA* trs = new AliITStrackSA();
572 fPoint1[0]=primaryVertex[0];
573 fPoint1[1]=primaryVertex[1];
574 fPoint2[0]=fTable->GetXCluster(1,ncl2);
575 fPoint2[1]=fTable->GetYCluster(1,ncl2);
577 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
578 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
579 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
580 trs,primaryVertex[2],pflag,fTable);
588 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
589 if(nn[nnp]!=0) fl+=1;
592 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert,opt);
595 Int_t nct = trs->GetNumberOfClustersSA();
597 Int_t index = trs->GetClusterIndexSA(nct);
598 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
599 if(kl->IsUsed()==1) kl->Use();
604 AliESDtrack outtrack;
605 outtrack.UpdateTrackParams(tr2,AliESDtrack::kITSin);
606 event->AddTrack(&outtrack);
608 Int_t nct = tr2->GetNumberOfClusters();
610 Int_t index = tr2->GetClusterIndex(nct);
611 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
613 if(kl->IsUsed()==0) kl->Use();
617 Int_t nct = trs->GetNumberOfClustersSA();
619 Int_t index = trs->GetClusterIndexSA(nct);
620 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
622 if(kl->IsUsed()==1) kl->Use();
626 }//end loop on clusters of layer2
634 Info("FindTracks","Number of found tracks: %d",event->GetNumberOfTracks());
639 //________________________________________________________________________
641 AliITStrackV2* AliITStrackerSA::FitTrack(AliITStrackSA* tr,Double_t *primaryVertex,Double_t *errorsprimvert,char *opt){
642 //fit of the found track
647 Bool_t sixpoints= choice.Contains("6/6");
649 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
650 for(Int_t i=0;i<fGeom->GetNlayers();i++){
651 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
653 TObjArray* listoftracks=new TObjArray(0,0);
654 AliITStrackV2* otrack2;
655 Int_t nclusters = tr->GetNumberOfClustersSA();
656 TObjArray** listlayer = new TObjArray*[fGeom->GetNlayers()];
657 for(Int_t i=0;i<fGeom->GetNlayers();i++){
658 listlayer[i] = new TObjArray(0,0);
668 Int_t * nnn = new Int_t[fGeom->GetNlayers()];
669 for(Int_t i=0;i<fGeom->GetNlayers();i++)nnn[i]=0;
671 for(Int_t ncl=0;ncl<nclusters;ncl++){
672 Int_t index = tr->GetClusterIndexSA(ncl);
673 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(index);
675 if(cl->IsUsed()==1) cl->Use();
676 Int_t lay = (index & 0xf0000000) >> 28;
677 if(lay==0) { listlayer[0]->AddLast(cl); clind0[nnn[0]]=index;nnn[0]++;}
678 if(lay==1) { listlayer[1]->AddLast(cl); clind1[nnn[1]]=index;nnn[1]++;}
679 if(lay==2) { listlayer[2]->AddLast(cl); clind2[nnn[2]]=index;nnn[2]++;}
680 if(lay==3) { listlayer[3]->AddLast(cl); clind3[nnn[3]]=index;nnn[3]++;}
681 if(lay==4) { listlayer[4]->AddLast(cl); clind4[nnn[4]]=index;nnn[4]++;}
682 if(lay==5) { listlayer[5]->AddLast(cl); clind5[nnn[5]]=index;nnn[5]++;}
687 Int_t * end = new Int_t[fGeom->GetNlayers()];
688 for(Int_t i=0;i<fGeom->GetNlayers();i++){
689 if(listlayer[i]->GetEntries()==0) end[i]=1;
690 else end[i]=listlayer[i]->GetEntries();
693 for(Int_t l1=0;l1<end[0];l1++){//loop on layer 1
694 AliITSclusterV2* cl0 = (AliITSclusterV2*)listlayer[0]->At(l1);
695 TVector3** recp = new TVector3*[3];
696 TVector3** errs = new TVector3*[3];
697 recp[0] = new TVector3(primaryVertex[0],primaryVertex[1],primaryVertex[2]);
698 errs[0] = new TVector3(errorsprimvert[0],errorsprimvert[1],errorsprimvert[2]);
699 Double_t x1,y1,z1,sx1,sy1,sz1;
700 Double_t x2,y2,z2,sx2,sy2,sz2;
701 AliITSclusterV2* p1=0;
702 AliITSclusterV2* p2=0;
703 Int_t index1=clind0[l1];
705 for(Int_t l2=0;l2<end[1];l2++){//loop on layer 2
706 AliITSclusterV2* cl1 = (AliITSclusterV2*)listlayer[1]->At(l2);
708 for(Int_t l3=0;l3<end[2];l3++){ //loop on layer 3
709 AliITSclusterV2* cl2 = (AliITSclusterV2*)listlayer[2]->At(l3);
711 if(cl0==0 && cl1!=0) {
712 p2 = cl2;index1=clind2[l3];
716 if(cl0!=0 && cl1==0){
718 p2=cl2;index2=clind2[l3];
720 if(cl0!=0 && cl1!=0){
724 Int_t lay1=(index1 & 0xf0000000) >> 28;
725 Int_t cln1=(index1 & 0x0fffffff) >> 00;
726 Int_t lay2=(index2 & 0xf0000000) >> 28;
727 Int_t cln2=(index2 & 0x0fffffff) >> 00;
728 x1 = fTable->GetXCluster(lay1,cln1);
729 x2 = fTable->GetXCluster(lay2,cln2);
730 y1 = fTable->GetYCluster(lay1,cln1);
731 y2 = fTable->GetYCluster(lay2,cln2);
732 z1 = fTable->GetZCluster(lay1,cln1);
733 z2 = fTable->GetZCluster(lay2,cln2);
734 sx1 = fTable->GetXClusterError(lay1,cln1);
735 sx2 = fTable->GetXClusterError(lay2,cln2);
736 sy1 = fTable->GetYClusterError(lay1,cln1);
737 sy2 = fTable->GetYClusterError(lay2,cln2);
738 sz1 = fTable->GetZClusterError(lay1,cln1);
739 sz2 = fTable->GetZClusterError(lay2,cln2);
740 Double_t phi1 = fTable->GetPhiCluster(lay1,cln1);
741 Int_t module1 = p1->GetDetectorIndex()+firstmod[0];
742 recp[1] = new TVector3(x1,y1,z1);
743 errs[1] = new TVector3(sx1,sy1,sz1);
744 recp[2] = new TVector3(x2,y2,z2);
745 errs[2] = new TVector3(sx2,sy2,sz2);
747 //fit on the Riemann sphere
748 Float_t seed1,seed2,seed3;
749 AliITSRiemannFit fit;
750 Int_t rf = fit.FitHelix(3,recp,errs,seed1,seed2,seed3); //this gives phi,tgl,curvature to start Kalman Filter
756 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
757 phi=seed1+1.5*TMath::Pi();
760 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
761 phi=seed1+0.5*TMath::Pi();
764 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
765 phi=seed1-0.5*TMath::Pi();
770 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
771 phi=seed1+0.5*TMath::Pi();
774 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
775 phi=seed1-0.5*TMath::Pi();
778 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
779 phi=seed1-1.5*TMath::Pi();
784 Int_t layer,ladder,detector;
785 fGeom->GetModuleId(module1,layer,ladder,detector);
786 Float_t yclu1 = p1->GetY();
787 Float_t zclu1 = p1->GetZ();
788 Double_t cv=Curvature(primaryVertex[0],primaryVertex[1],x1,y1,x2,y2);
790 for(Int_t l4=0;l4<end[3];l4++){ //loop on layer 4
791 AliITSclusterV2* cl3 = (AliITSclusterV2*)listlayer[3]->At(l4);
792 for(Int_t l5=0;l5<end[4];l5++){ //loop on layer 5
793 AliITSclusterV2* cl4 = (AliITSclusterV2*)listlayer[4]->At(l5);
794 for(Int_t l6=0;l6<end[5];l6++){ //loop on layer 6
795 AliITSclusterV2* cl5 = (AliITSclusterV2*)listlayer[5]->At(l6);
796 AliITStrackSA* trac = new AliITStrackSA(layer,ladder,detector,yclu1,zclu1,phi,tgl,cv,1);
798 if(cl5!=0) trac->AddClusterV2(5,(clind5[l6] & 0x0fffffff)>>0);
799 if(cl4!=0) trac->AddClusterV2(4,(clind4[l5] & 0x0fffffff)>>0);
800 if(cl3!=0) trac->AddClusterV2(3,(clind3[l4] & 0x0fffffff)>>0);
801 if(cl2!=0) trac->AddClusterV2(2,(clind2[l3] & 0x0fffffff)>>0);
802 if(cl1!=0) trac->AddClusterV2(1,(clind1[l2] & 0x0fffffff)>>0);
803 if(cl0!=0) trac->AddClusterV2(0,(clind0[l1] & 0x0fffffff)>>0);
805 //fit with Kalman filter using AliITStrackerV2::RefitAt()
807 AliITStrackV2* ot = new AliITStrackV2(*trac);
809 ot->ResetCovariance();
812 if(RefitAt(49.,ot,trac)){ //fit from layer 1 to layer 6
814 otrack2 = new AliITStrackV2(*ot);
815 otrack2->ResetCovariance();
816 otrack2->ResetClusters();
817 //fit from layer 6 to layer 1
818 if(RefitAt(3.7,otrack2,ot)) listoftracks->AddLast(otrack2);
828 for(Int_t i=1;i<3;i++){
842 Int_t dim=listoftracks->GetEntries();
845 for(Int_t i=0;i<fGeom->GetNlayers();i++){
852 AliITStrackV2* otrack =(AliITStrackV2*)FindTrackLowChiSquare(listoftracks,dim);
854 if(otrack==0) return 0;
855 Int_t * indexc = new Int_t[fGeom->GetNlayers()];
856 for(Int_t i=0;i<fGeom->GetNlayers();i++) indexc[i]=0;
857 for(Int_t nind=0;nind<otrack->GetNumberOfClusters();nind++){
858 indexc[nind] = otrack->GetClusterIndex(nind);
860 AliITSclusterV2* cl0 = (AliITSclusterV2*)GetCluster(indexc[0]);
861 AliITSclusterV2* cl1 = (AliITSclusterV2*)GetCluster(indexc[1]);
862 AliITSclusterV2* cl2 = (AliITSclusterV2*)GetCluster(indexc[2]);
863 AliITSclusterV2* cl3 = (AliITSclusterV2*)GetCluster(indexc[3]);
864 AliITSclusterV2* cl4 = (AliITSclusterV2*)GetCluster(indexc[4]);
865 Int_t labl[3]={-1,-1,-1};
866 if(otrack->GetNumberOfClusters()==fGeom->GetNlayers()){
867 AliITSclusterV2* cl5 = (AliITSclusterV2*)GetCluster(indexc[5]);
868 labl[0]=cl5->GetLabel(0);
869 labl[1]=cl5->GetLabel(1);
870 labl[2]=cl5->GetLabel(2);
873 if(otrack->GetNumberOfClusters()==(fGeom->GetNlayers()-1)){
878 Int_t numberofpoints;
879 if(sixpoints) numberofpoints=6;
880 else numberofpoints=5;
881 Int_t label = Label(cl0->GetLabel(0),cl1->GetLabel(0),
882 cl2->GetLabel(0),cl3->GetLabel(0),
883 cl4->GetLabel(0),labl[0],
884 cl0->GetLabel(1),cl1->GetLabel(1),
885 cl2->GetLabel(1),cl3->GetLabel(1),
886 cl4->GetLabel(1),labl[1],
887 cl0->GetLabel(2),cl1->GetLabel(2),
888 cl2->GetLabel(2),cl3->GetLabel(2),
889 cl4->GetLabel(2),labl[2],numberofpoints);
891 otrack->SetLabel(label);
893 for(Int_t i=0;i<fGeom->GetNlayers();i++){
902 //_______________________________________________________________________
903 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,TTree* treev2, TTree* clustertree){
904 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
905 //(or AliITStrackV2 tracks found with function FindTracks of this class)
909 if(fVert)delete fVert;
910 fVert = fVertexer->FindVertexForCurrentEvent(evnum);
913 gAlice->GetEvent(evnum);
915 Fatal("FindTracks","Vertex is missing\n");
919 Double_t primaryVertex[3];
920 fVert->GetXYZ(primaryVertex);
923 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
924 fTable->FillArray(clustertree,evnum);
925 fTable->FillArrayCoorAngles(clustertree,evnum);
927 SetEventNumber(evnum);
928 if(GetFlagLoadedClusters()==0) LoadClusters(clustertree);
929 SetFlagLoadedClusters(1);
931 TBranch* bra = (TBranch*)treev2->GetBranch("tracks");
932 if(!bra) Warning("UseFoundTracksV2","No branch for track tree");
933 AliITStrackV2* ttrrt = new AliITStrackV2;
934 bra->SetAddress(&ttrrt);
936 for(Int_t nj=0;nj<treev2->GetEntries();nj++){
937 treev2->GetEvent(nj);
938 Int_t ncl = ttrrt->GetNumberOfClusters();
939 for(Int_t k=0;k<ncl;k++){
940 Int_t index = ttrrt->GetClusterIndex(k);
941 AliITSclusterV2* clui = (AliITSclusterV2*)GetCluster(index);
942 if(clui->IsUsed()==0) clui->Use();
949 //_______________________________________________________________________
950 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,AliESD *event, TTree*
952 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
956 if(fVert)delete fVert;
957 fVert = fVertexer->FindVertexForCurrentEvent(evnum);
960 gAlice->GetEvent(evnum);
962 Fatal("FindTracks","Vertex is missing\n");
966 Double_t primaryVertex[3];
967 fVert->GetXYZ(primaryVertex);
970 fTable = new AliITSclusterTable(fGeom,this,primaryVertex);
971 fTable->FillArray(clustertree,evnum);
972 fTable->FillArrayCoorAngles(clustertree,evnum);
974 SetEventNumber(evnum);
975 if(GetFlagLoadedClusters()==0) LoadClusters(clustertree);
976 SetFlagLoadedClusters(1);
977 Int_t ntracks = event->GetNumberOfTracks();
979 AliESDtrack *esd=event->GetTrack(ntracks);
980 if ((esd->GetStatus()&
981 AliESDtrack::kITSin|AliESDtrack::kTPCin)==0) continue;
983 Int_t ncl = esd->GetITSclusters(idx);
984 for(Int_t clu=0; clu<ncl; clu++){
985 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(idx[clu]);
986 if(cl->IsUsed()==0) cl->Use();
990 Info("UseFoundTracksV2","Clusters of tracks prolonged from TPC deleted");
995 //_______________________________________________________
996 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Int_t pflag,AliITSclusterTable* table){
997 //function used to to find the clusters associated to the track
999 AliITSlayer &lay = fgLayers[layer];
1000 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
1001 for(Int_t i=0;i<fGeom->GetNlayers();i++){
1002 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
1006 Float_t cx1,cx2,cy1,cy2;
1007 FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
1008 Int_t fun = FindIntersection(fCoef1,fCoef2,fCoef3,-(lay.GetR()*lay.GetR()),cx1,cy1,cx2,cy2);
1009 if(fun==0) return 0;
1010 Double_t fi1 =TMath::ATan2(cy1,cx1);
1011 Double_t fi2 =TMath::ATan2(cy2,cx2);
1012 fPhiEstimate = ChoosePoint(fi1,fi2,fPhic);
1015 Double_t zed = TMath::Tan(fLambdac)*lay.GetR()+zvertex;
1016 Double_t zed1 = TMath::Tan(fLambdac+lambdawindow)*lay.GetR()+zvertex;
1017 Double_t zed2 = TMath::Tan(fLambdac-lambdawindow)*lay.GetR()+zvertex;
1019 Double_t fi = fPhiEstimate;
1020 Int_t nmod = lay.FindDetectorIndex(fi,zed)+firstmod[layer];
1022 Int_t nm[8]={0,0,0,0,0,0,0,0};
1023 nm[0] = lay.FindDetectorIndex(fi+phiwindow,zed)+firstmod[layer];
1024 nm[1] = lay.FindDetectorIndex(fi-phiwindow,zed)+firstmod[layer];
1025 nm[2] = lay.FindDetectorIndex(fi,zed1)+firstmod[layer];
1026 nm[3] = lay.FindDetectorIndex(fi,zed2)+firstmod[layer];
1027 nm[4] = lay.FindDetectorIndex(fi+phiwindow,zed1)+firstmod[layer];
1028 nm[5] = lay.FindDetectorIndex(fi-phiwindow,zed1)+firstmod[layer];
1029 nm[6] = lay.FindDetectorIndex(fi+phiwindow,zed2)+firstmod[layer];
1030 nm[7] = lay.FindDetectorIndex(fi-phiwindow,zed2)+firstmod[layer];
1033 if(nmod<0) return 0;
1035 TArrayI* array =(TArrayI*)table->GetListOfClusters(nmod);
1036 TArrayI* list = new TArrayI(array->GetSize());
1037 for(Int_t i=0;i<array->GetSize();i++){
1038 Int_t in=(Int_t)array->At(i);
1043 for(Int_t ii=0;ii<8;ii++){
1044 if(nm[ii]!=nmod && nm[ii]>=0){
1045 TArrayI* ar =(TArrayI*)table->GetListOfClusters(nm[ii]);
1046 list->Set(list->GetSize()+ar->GetSize());
1047 for(Int_t j=0;j<ar->GetSize();j++){
1048 Int_t in=(Int_t)ar->At(j);
1054 for(Int_t i=0;i<list->GetSize();i++){
1055 Int_t index = (Int_t)list->At(i);
1056 AliITSclusterV2* cllay = lay.GetCluster(index);
1057 if(cllay==0) continue;
1058 if(cllay->IsUsed()==1) continue;
1059 Double_t phi = fTable->GetPhiCluster(layer,index);
1060 Double_t lambda= fTable->GetLambdaCluster(layer,index);
1062 if(TMath::Abs(fLambdac-lambda)<lambdawindow &&
1063 TMath::Abs(fPhiEstimate-phi)<phiwindow){
1066 if(trs->GetNumberOfClustersSA()==20){
1070 trs->AddClusterSA(layer,index);
1073 fPointc[0]=fTable->GetXCluster(layer,index);
1074 fPointc[1]=fTable->GetYCluster(layer,index);
1086 //________________________________________________________________
1087 void AliITStrackerSA::UpdatePoints(){
1088 //update of points for the estimation of the curvature
1090 fPoint1[0]=fPoint2[0];
1091 fPoint1[1]=fPoint2[1];
1092 fPoint2[0]=fPoint3[0];
1093 fPoint2[1]=fPoint3[1];
1094 fPoint3[0]=fPointc[0];
1095 fPoint3[1]=fPointc[1];
1102 //___________________________________________________________________
1103 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){
1105 //given (x,y) of three recpoints (in global coordinates)
1106 //returns the parameters a,b,c of circonference x*x + y*y +a*x + b*y +c
1108 Float_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1109 if(den==0) return 0;
1110 a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1111 b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1112 c = -x1*x1-y1*y1-a*x1-b*y1;
1115 //__________________________________________________________________________
1116 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){
1118 //Finds the intersection between the circonference of the track and the circonference centered in (0,0) represented by one layer
1119 //c2 is -rlayer*rlayer
1123 Float_t aA = (b1*b1)/(a1*a1)+1;
1124 Float_t bB = (-2*m*b1/(a1*a1));
1125 Float_t cC = c2+(m*m)/(a1*a1);
1126 if((bB*bB-4*aA*cC)<0) return 0;
1128 y1 = (-bB+TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA);
1129 y2 = (-bB-TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA);
1130 x1 = (c2-c1-b1*y1)/a1;
1131 x2 = (c2-c1-b1*y2)/a1;
1135 //____________________________________________________________________
1136 Double_t AliITStrackerSA::Curvature(Double_t x1,Double_t y1,Double_t
1137 x2,Double_t y2,Double_t x3,Double_t y3){
1139 //calculates the curvature of track
1140 Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
1141 if(den==0) return 0;
1142 Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
1143 Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
1144 Double_t c = -x1*x1-y1*y1-a*x1-b*y1;
1147 if((a*a+b*b-4*c)<0) return 0;
1148 Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.;
1149 if(rad==0) return 0;
1151 if((x1>0 && y1>0 && x1<xc)) rad*=-1;
1152 if((x1<0 && y1>0 && x1<xc)) rad*=-1;
1153 // if((x1<0 && y1<0 && x1<xc)) rad*=-1;
1154 // if((x1>0 && y1<0 && x1<xc)) rad*=-1;
1159 //____________________________________________________________________
1160 Double_t AliITStrackerSA::ChoosePoint(Double_t p1, Double_t p2, Double_t pp){
1162 //Returns the point closest to pp
1164 Double_t diff1 = p1-pp;
1165 Double_t diff2 = p2-pp;
1167 if(TMath::Abs(diff1)<TMath::Abs(diff2)) fPhiEstimate=p1;
1168 else fPhiEstimate=p2;
1169 return fPhiEstimate;
1174 //_________________________________________________________________
1175 AliITStrackV2* AliITStrackerSA::FindTrackLowChiSquare(TObjArray* tracklist, Int_t dim) const {
1176 // returns track with lowes chi square
1178 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(0);
1181 if(dim==0) return 0;
1182 Double_t * chi2 = new Double_t[dim];
1183 Int_t * index = new Int_t[dim];
1184 for(Int_t i=0;i<dim;i++){
1185 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(i);
1186 chi2[i]=trk->GetChi2();
1190 Int_t w=0;Double_t value;
1193 for(Int_t j=w+1;j<dim;j++){
1194 if(chi2[w]<chi2[j]){
1206 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(index[dim-1]);
1213 //__________________________________________________________
1214 Int_t AliITStrackerSA::FindLabel(Int_t l1, Int_t l2, Int_t l3, Int_t l4, Int_t l5, Int_t l6){
1216 //function used to determine the track label
1218 Int_t lb[6] = {l1,l2,l3,l4,l5,l6};
1219 Int_t aa[6]={1,1,1,1,1,1};
1222 Int_t k=0;Int_t w=0;Int_t num=6;
1223 if(lb[5]==-1) num=5;
1227 for(Int_t i=k+1;i<num;i++){
1229 if(lb[k]==lb[i] && aa[k]!=0){
1240 for(Int_t j=0;j<6;j++){
1252 if(num==6) return lb[5];
1256 //_____________________________________________________________________________
1257 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,
1258 Int_t gl12, Int_t gl13, Int_t gl14,Int_t gl15, Int_t gl16, Int_t gl17, Int_t gl18, Int_t numberofpoints){
1261 //function used to assign label to the found track. If track is fake, the label is negative
1263 Int_t lb0[6] = {gl1,gl2,gl3,gl4,gl5,gl6};
1264 Int_t lb1[6] = {gl7,gl8,gl9,gl10,gl11,gl12};
1265 Int_t lb2[6] = {gl13,gl14,gl15,gl16,gl17,gl18};
1266 Int_t ll=FindLabel(lb0[0],lb0[1],lb0[2],lb0[3],lb0[4],lb0[5]);
1267 Int_t lflag=0;Int_t num=6;
1268 if(lb0[5]==-1 && lb1[5]==-1 && lb2[5]==-1) num=5;
1270 for(Int_t i=0;i<num;i++){
1271 if(lb0[i]==ll || lb1[i]==ll || lb2[i]==ll) lflag+=1;
1274 if(lflag>=numberofpoints) return ll;
1280 //_____________________________________________________________________________
1281 void AliITStrackerSA::SetWindowSizes(Int_t n, Double_t *phi, Double_t *lam){
1282 // Set sizes of the phi and lambda windows used for track finding
1284 if(phi){ // user defined values
1285 fPhiWin = new Double_t[fNloop];
1286 fLambdaWin = new Double_t[fNloop];
1287 for(Int_t k=0;k<fNloop;k++){
1289 fLambdaWin[k]=lam[k];
1292 else { // default values
1294 Double_t phid[40] = {0.001,0.0015,0.002,0.0023,0.0025,0.0027,0.003,
1295 0.0033,0.0035,0.0037,0.004,0.0043,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.0105,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1301 Double_t lambdad[40] = {0.001,0.0015,0.002,0.0023,0.0025,0.0027,0.003,
1302 0.0033,0.0035,0.0037,0.004,0.0043,0.0045,0.0047,
1303 0.005,0.0053,0.0055,
1304 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1305 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1306 0.01,0.015,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1311 Warning("SetWindowSizes","Number of loop forced to the default value %d",fNloop);
1314 fPhiWin = new Double_t[fNloop];
1315 fLambdaWin = new Double_t[fNloop];
1317 for(Int_t k=0;k<fNloop;k++){
1319 fLambdaWin[k]=lambdad[k];