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"
40 ClassImp(AliITStrackerSA)
42 //____________________________________________________________________________
43 AliITStrackerSA::AliITStrackerSA():AliITStrackerV2(){
44 // Default constructor
47 //____________________________________________________________________________
48 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertex *vert):AliITStrackerV2(geom)
50 // Standard constructor (Vertex is known and passed to this obj.)
56 //______________________________________________________________________
57 AliITStrackerSA::AliITStrackerSA(const AliITStrackerSA &trkr) :
58 AliITStrackerV2(trkr) {
60 // Copies are not allowed. The method is protected to avoid misuse.
61 Error("AliITStrackerSA","Copy constructor not allowed\n");
64 //______________________________________________________________________
65 AliITStrackerSA& AliITStrackerSA::operator=(const
66 AliITStrackerSA& /* trkr */){
67 // Assignment operator
68 // Assignment is not allowed. The method is protected to avoid misuse.
69 Error("= operator","Assignment operator not allowed\n");
73 //____________________________________________________________________________
74 AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertexer *vertexer):AliITStrackerV2(geom)
76 // Standard constructor (Vertex is unknown - vertexer is passed to this obj)
82 //____________________________________________________________________________
83 AliITStrackerSA::AliITStrackerSA(AliITStrackerSA& tracker):AliITStrackerV2(){
85 fPhiEstimate = tracker.fPhiEstimate;
86 for(Int_t i=0;i<2;i++){
87 fPoint1[i]=tracker.fPoint1[i];
88 fPoint2[i]=tracker.fPoint2[i];
89 fPoint3[i]=tracker.fPoint3[i];
90 fPointc[i]=tracker.fPointc[i];
92 fLambdac = tracker.fLambdac;
93 fPhic = tracker.fPhic;
94 fCoef1 = tracker.fCoef1;
95 fCoef2 = tracker.fCoef2;
96 fCoef3 = tracker.fCoef3;
97 fNloop = tracker.fNloop;
98 fPhiWin = tracker.fPhiWin;
99 fLambdaWin = tracker.fLambdaWin;
100 if(tracker.fVertexer && tracker.fVert){
101 fVert = new AliITSVertex(*tracker.fVert);
104 fVert = tracker.fVert;
106 fVertexer = tracker.fVertexer;
107 fGeom = tracker.fGeom;
108 fFlagLoad = tracker.fFlagLoad;
109 fTable = tracker.fTable;
112 //____________________________________________________________________________
113 AliITStrackerSA::~AliITStrackerSA(){
115 // if fVertexer is not null, the AliITSVertex obj. is owned by this class
116 // and is deleted here
118 if(fVert)delete fVert;
124 if(fPhiWin)delete []fPhiWin;
125 if(fLambdaWin)delete []fLambdaWin;
129 //____________________________________________________________________________
130 void AliITStrackerSA::Init(){
131 // Reset all data members
133 for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
148 //_______________________________________________________________________
149 void AliITStrackerSA::ResetForFinding(){
150 // Reset data members used in all loops during track finding
152 for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;}
161 //____________________________________________________________________________
162 void AliITStrackerSA::FindTracks(TTree *clusterTree, TTree *out,Int_t evnumber,char *opt){
164 /**************************************************************************
165 * Options: to use only ITS execute FindTracks *
166 * to define good tracks with 6 points out of 6 in the ITS write *
167 * *opt="6/6". The default is *opt="6/6" *
170 * Example: to execute function with only the ITS (no combined tracking *
171 * with TPC+ITS) and requiring 5/6 points to define a good track *
172 * use: FindTracks(treein,treeout,evnumber,"5/6") *
173 * to execute combined tracking, before using FindTracks use *
175 *************************************************************************/
179 Bool_t sixpoints= choice.Contains("6/6");
181 //Fill array with cluster indices for each module
183 fTable = new AliITSclusterTable(fGeom,this);
184 fTable->FillArray(clusterTree,evnumber);
188 if(fVert)delete fVert;
189 fVert = fVertexer->FindVertexForCurrentEvent(evnumber);
192 gAlice->GetEvent(evnumber);
194 Fatal("FindTracks","Vertex is missing\n");
198 Double_t primaryVertex[3];
199 Double_t errorsprimvert[3];
200 fVert->GetXYZ(primaryVertex);
201 fVert->GetSigmaXYZ(errorsprimvert);
202 if(errorsprimvert[0]==0 || errorsprimvert[1]==0){
203 Warning("FindTracks","Set errors on vertex positions x and y at 0.0001");
204 errorsprimvert[0]=0.0001;
205 errorsprimvert[1]=0.0001;
207 fVert->PrintStatus();
209 SetEventNumber(evnumber);
210 if(GetFlagLoadedClusters()==0) LoadClusters(clusterTree);
212 //Fill tree for found tracks
213 AliITStrackV2* outrack=0;
214 TBranch* branch=out->Branch("Tracks","AliITStrackV2",&outrack,32000,0);
215 if (!branch) out->Branch("tracks","AliITStrackV2",&outrack,32000,3);
216 else branch->SetAddress(&outrack);
219 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
220 for(Int_t i=0;i<fGeom->GetNlayers();i++){
221 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
223 // firstmod [i] number of the first module in the ITS layer i.
225 AliITSlayer &layer=fgLayers[0]; // first layer
227 //loop on the different windows
228 for(Int_t nloop=0;nloop<fNloop;nloop++){
229 Int_t ncl=layer.GetNumberOfClusters();
230 while(ncl--){ //loop starting from layer 0
233 AliITSclusterV2* cl = layer.GetCluster(ncl);
234 if(cl->IsUsed()==1) continue;
236 Int_t module = cl->GetDetectorIndex()+firstmod[0];
237 GetClusterGCoordinates(cl,module,x,y,z);
238 fPhic = TMath::ATan2(y,x);
239 fLambdac = TMath::ATan2(z-primaryVertex[2],TMath::Sqrt((x-primaryVertex[0])*(x-primaryVertex[0])+(y-primaryVertex[1])*(y-primaryVertex[1])));
240 fPhiEstimate = fPhic;
241 AliITStrackSA* trs = new AliITStrackSA();
242 fPoint1[0]=primaryVertex[0];
243 fPoint1[1]=primaryVertex[1];
247 Int_t * nn = new Int_t[fGeom->GetNlayers()];//counter for clusters on each layer
248 for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;}
249 nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable);
250 nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable);
254 fPoint3[0] = fPointc[0];
255 fPoint3[1] = fPointc[1];
257 nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable);
258 if(nn[1]==0 && nn[2]==0) pflag=0;
259 if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();}
260 if(nn[2]!=0 && nn[1]==0){
262 fPoint3[0]=fPointc[0];
263 fPoint3[1]=fPointc[1];
266 nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable);
268 if(nn[3]!=0) UpdatePoints();
269 nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable);
271 if(nn[4]!=0) UpdatePoints();
272 nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable);
276 Int_t numberofpoints;
277 if(sixpoints) numberofpoints=6; //check of the candidate track
278 else numberofpoints=5; //if track is good (with the required number
279 for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file
280 if(nn[nnp]!=0) layOK+=1;
282 if(layOK>=numberofpoints){
283 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert,opt);
285 Int_t nct = trs->GetNumberOfClustersSA();
287 Int_t index = trs->GetClusterIndexSA(nct);
288 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
289 if(kl->IsUsed()==1) kl->Use();
296 Int_t nct = tr2->GetNumberOfClusters();
299 Int_t index = tr2->GetClusterIndex(nct);
300 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
301 if(kl->IsUsed()==0) kl->Use();
306 Int_t nct = trs->GetNumberOfClustersSA();
308 Int_t index = trs->GetClusterIndexSA(nct);
309 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
310 if(kl->IsUsed()==1) kl->Use();
315 }//end loop on clusters of layer1
319 //if 5/6 points are required, second loop starting
320 //from second layer, to find tracks with point of
324 // counter for clusters on each layer
325 Int_t * nn = new Int_t[fGeom->GetNlayers()-1];
326 for(Int_t nloop=0;nloop<fNloop;nloop++){
327 AliITSlayer &layer2=fgLayers[1]; //loop on layer 2
328 Int_t ncl2=layer2.GetNumberOfClusters();
329 while(ncl2--){ //loop starting from layer 2
332 AliITSclusterV2* cl = layer2.GetCluster(ncl2);
333 if(cl->IsUsed()==1) continue;
335 Int_t module = cl->GetDetectorIndex()+firstmod[1];
336 GetClusterGCoordinates(cl,module,x,y,z);
337 fPhic = TMath::ATan2(y,x);
338 fLambdac = TMath::ATan2(z-primaryVertex[2],TMath::Sqrt((x-primaryVertex[0])*(x-primaryVertex[0])+(y-primaryVertex[1])*(y-primaryVertex[1])));
340 fPhiEstimate = fPhic;
341 AliITStrackSA* trs = new AliITStrackSA();
343 fPoint1[0]=primaryVertex[0];
344 fPoint1[1]=primaryVertex[1];
348 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0;
349 for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){
350 nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop],
351 trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable);
359 for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){
360 if(nn[nnp]!=0) fl+=1;
363 AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert,opt);
365 Int_t nct = trs->GetNumberOfClustersSA();
367 Int_t index = trs->GetClusterIndexSA(nct);
368 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
369 if(kl->IsUsed()==1) kl->Use();
375 Int_t nct = tr2->GetNumberOfClusters();
377 Int_t index = tr2->GetClusterIndex(nct);
378 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
380 if(kl->IsUsed()==0) kl->Use();
384 Int_t nct = trs->GetNumberOfClustersSA();
386 Int_t index = trs->GetClusterIndexSA(nct);
387 AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index);
389 if(kl->IsUsed()==1) kl->Use();
393 }//end loop on clusters of layer2
402 //_______________________________________________________________________
403 void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,TTree* treev2, TTree* clustertree){
404 // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking
406 fTable = new AliITSclusterTable(fGeom,this);
407 fTable->FillArray(clustertree,evnum);
409 LoadClusters(clustertree);
410 SetFlagLoadedClusters(1);
412 TBranch* bra = (TBranch*)treev2->GetBranch("tracks");
413 if(!bra) Warning("UseFoundTracksV2","No branch for track tree");
414 AliITStrackV2* ttrrt = new AliITStrackV2;
415 bra->SetAddress(&ttrrt);
417 for(Int_t nj=0;nj<treev2->GetEntries();nj++){
418 treev2->GetEvent(nj);
419 Int_t ncl = ttrrt->GetNumberOfClusters();
420 for(Int_t k=0;k<ncl;k++){
421 Int_t index = ttrrt->GetClusterIndex(k);
422 AliITSclusterV2* clui = (AliITSclusterV2*)GetCluster(index);
423 if(clui->IsUsed()==0) clui->Use();
430 //________________________________________________________________________
432 AliITStrackV2* AliITStrackerSA::FitTrack(AliITStrackSA* tr,Double_t *primaryVertex,Double_t *errorsprimvert,char *opt){
433 //fit of the found track
437 Bool_t sixpoints= choice.Contains("6/6");
439 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
440 for(Int_t i=0;i<fGeom->GetNlayers();i++){
441 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
443 TObjArray* listoftracks=new TObjArray(0,0);
444 AliITStrackV2* otrack2;
445 Int_t nclusters = tr->GetNumberOfClustersSA();
446 TObjArray** listlayer = new TObjArray*[fGeom->GetNlayers()];
447 for(Int_t i=0;i<fGeom->GetNlayers();i++){
448 listlayer[i] = new TObjArray(0,0);
458 Int_t * nnn = new Int_t[fGeom->GetNlayers()];
459 for(Int_t i=0;i<fGeom->GetNlayers();i++)nnn[i]=0;
461 for(Int_t ncl=0;ncl<nclusters;ncl++){
462 Int_t index = tr->GetClusterIndexSA(ncl);
463 AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(index);
464 if(cl->IsUsed()==1) cl->Use();
465 Int_t lay = (index & 0xf0000000) >> 28;
466 if(lay==0) { listlayer[0]->AddLast(cl); clind0[nnn[0]]=index;nnn[0]++;}
467 if(lay==1) { listlayer[1]->AddLast(cl); clind1[nnn[1]]=index;nnn[1]++;}
468 if(lay==2) { listlayer[2]->AddLast(cl); clind2[nnn[2]]=index;nnn[2]++;}
469 if(lay==3) { listlayer[3]->AddLast(cl); clind3[nnn[3]]=index;nnn[3]++;}
470 if(lay==4) { listlayer[4]->AddLast(cl); clind4[nnn[4]]=index;nnn[4]++;}
471 if(lay==5) { listlayer[5]->AddLast(cl); clind5[nnn[5]]=index;nnn[5]++;}
476 Int_t * end = new Int_t[fGeom->GetNlayers()];
477 for(Int_t i=0;i<fGeom->GetNlayers();i++){
478 if(listlayer[i]->GetEntries()==0) end[i]=1;
479 else end[i]=listlayer[i]->GetEntries();
482 for(Int_t l1=0;l1<end[0];l1++){//loop on layer 1
483 AliITSclusterV2* cl0 = (AliITSclusterV2*)listlayer[0]->At(l1);
484 TVector3** recp = new TVector3*[3];
485 TVector3** errs = new TVector3*[3];
486 recp[0] = new TVector3(primaryVertex[0],primaryVertex[1],primaryVertex[2]);
487 errs[0] = new TVector3(errorsprimvert[0],errorsprimvert[1],errorsprimvert[2]);
488 Double_t x1,y1,z1,sx1,sy1,sz1;
489 Double_t x2,y2,z2,sx2,sy2,sz2;
490 AliITSclusterV2* p1=0;
491 AliITSclusterV2* p2=0;
493 for(Int_t l2=0;l2<end[1];l2++){//loop on layer 2
494 AliITSclusterV2* cl1 = (AliITSclusterV2*)listlayer[1]->At(l2);
495 for(Int_t l3=0;l3<end[2];l3++){ //loop on layer 3
496 AliITSclusterV2* cl2 = (AliITSclusterV2*)listlayer[2]->At(l3);
498 if(cl0==0 && cl1!=0) {
503 if(cl0!=0 && cl1==0){
507 if(cl0!=0 && cl1!=0){
511 Int_t module1 = p1->GetDetectorIndex()+firstmod[0];
512 Int_t module2 = p2->GetDetectorIndex()+firstmod[1];
513 GetClusterGCoordinates(p1,module1,x1,y1,z1);
514 GetClusterGCoordinates(p2,module2,x2,y2,z2);
515 GetClusterGErrors(p1,module1,sx1,sy1,sz1);
516 GetClusterGErrors(p2,module2,sx2,sy2,sz2);
517 Double_t phi1=TMath::ATan2(y1,x1);
518 recp[1] = new TVector3(x1,y1,z1);
519 errs[1] = new TVector3(sx1,sy1,sz1);
520 recp[2] = new TVector3(x2,y2,z2);
521 errs[2] = new TVector3(sx2,sy2,sz2);
523 //fit on the Riemann sphere
524 Float_t seed1,seed2,seed3;
525 AliITSRiemannFit fit;
526 Int_t rf = fit.FitHelix(3,recp,errs,seed1,seed2,seed3); //this gives phi,tgl,curvature to start Kalman Filter
532 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
533 phi=seed1+1.5*TMath::Pi();
536 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
537 phi=seed1+0.5*TMath::Pi();
540 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
541 phi=seed1-0.5*TMath::Pi();
546 if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){
547 phi=seed1+0.5*TMath::Pi();
550 if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){
551 phi=seed1-0.5*TMath::Pi();
554 if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){
555 phi=seed1-1.5*TMath::Pi();
560 Int_t layer,ladder,detector;
561 fGeom->GetModuleId(module1,layer,ladder,detector);
562 Float_t yclu1 = p1->GetY();
563 Float_t zclu1 = p1->GetZ();
564 Double_t cv=Curvature(primaryVertex[0],primaryVertex[1],x1,y1,x2,y2);
566 for(Int_t l4=0;l4<end[3];l4++){ //loop on layer 4
567 AliITSclusterV2* cl3 = (AliITSclusterV2*)listlayer[3]->At(l4);
568 for(Int_t l5=0;l5<end[4];l5++){ //loop on layer 5
569 AliITSclusterV2* cl4 = (AliITSclusterV2*)listlayer[4]->At(l5);
570 for(Int_t l6=0;l6<end[5];l6++){ //loop on layer 6
571 AliITSclusterV2* cl5 = (AliITSclusterV2*)listlayer[5]->At(l6);
572 AliITStrackSA* trac = new AliITStrackSA(layer,ladder,detector,yclu1,zclu1,phi,tgl,cv,1);
574 if(cl5!=0) trac->AddClusterV2(5,(clind5[l6] & 0x0fffffff)>>0);
575 if(cl4!=0) trac->AddClusterV2(4,(clind4[l5] & 0x0fffffff)>>0);
576 if(cl3!=0) trac->AddClusterV2(3,(clind3[l4] & 0x0fffffff)>>0);
577 if(cl2!=0) trac->AddClusterV2(2,(clind2[l3] & 0x0fffffff)>>0);
578 if(cl1!=0) trac->AddClusterV2(1,(clind1[l2] & 0x0fffffff)>>0);
579 if(cl0!=0) trac->AddClusterV2(0,(clind0[l1] & 0x0fffffff)>>0);
581 //fit with Kalman filter using AliITStrackerV2::RefitAt()
583 AliITStrackV2* ot = new AliITStrackV2(*trac);
585 ot->ResetCovariance();
588 if(RefitAt(49.,ot,trac)){ //fit from layer 1 to layer 6
590 otrack2 = new AliITStrackV2(*ot);
591 otrack2->ResetCovariance();
592 otrack2->ResetClusters();
593 //fit from layer 6 to layer 1
594 if(RefitAt(3.7,otrack2,ot)) listoftracks->AddLast(otrack2);
604 for(Int_t i=1;i<3;i++){
618 Int_t dim=listoftracks->GetEntries();
621 for(Int_t i=0;i<fGeom->GetNlayers();i++){
628 AliITStrackV2* otrack =(AliITStrackV2*)FindTrackLowChiSquare(listoftracks,dim);
630 if(otrack==0) return 0;
631 Int_t * indexc = new Int_t[fGeom->GetNlayers()];
632 for(Int_t i=0;i<fGeom->GetNlayers();i++) indexc[i]=0;
633 for(Int_t nind=0;nind<otrack->GetNumberOfClusters();nind++){
634 indexc[nind] = otrack->GetClusterIndex(nind);
636 AliITSclusterV2* cl0 = (AliITSclusterV2*)GetCluster(indexc[0]);
637 AliITSclusterV2* cl1 = (AliITSclusterV2*)GetCluster(indexc[1]);
638 AliITSclusterV2* cl2 = (AliITSclusterV2*)GetCluster(indexc[2]);
639 AliITSclusterV2* cl3 = (AliITSclusterV2*)GetCluster(indexc[3]);
640 AliITSclusterV2* cl4 = (AliITSclusterV2*)GetCluster(indexc[4]);
641 Int_t labl[3]={-1,-1,-1};
642 if(otrack->GetNumberOfClusters()==fGeom->GetNlayers()){
643 AliITSclusterV2* cl5 = (AliITSclusterV2*)GetCluster(indexc[5]);
644 labl[0]=cl5->GetLabel(0);
645 labl[1]=cl5->GetLabel(1);
646 labl[2]=cl5->GetLabel(2);
649 if(otrack->GetNumberOfClusters()==(fGeom->GetNlayers()-1)){
654 Int_t numberofpoints;
655 if(sixpoints) numberofpoints=6;
656 else numberofpoints=5;
657 Int_t label = Label(cl0->GetLabel(0),cl1->GetLabel(0),
658 cl2->GetLabel(0),cl3->GetLabel(0),
659 cl4->GetLabel(0),labl[0],
660 cl0->GetLabel(1),cl1->GetLabel(1),
661 cl2->GetLabel(1),cl3->GetLabel(1),
662 cl4->GetLabel(1),labl[1],
663 cl0->GetLabel(2),cl1->GetLabel(2),
664 cl2->GetLabel(2),cl3->GetLabel(2),
665 cl4->GetLabel(2),labl[2],numberofpoints);
667 otrack->SetLabel(label);
669 for(Int_t i=0;i<fGeom->GetNlayers();i++){
678 //________________________________________________________________
679 void AliITStrackerSA::UpdatePoints(){
680 //update of points for the estimation of the curvature
682 fPoint1[0]=fPoint2[0];
683 fPoint1[1]=fPoint2[1];
684 fPoint2[0]=fPoint3[0];
685 fPoint2[1]=fPoint3[1];
686 fPoint3[0]=fPointc[0];
687 fPoint3[1]=fPointc[1];
692 //_________________________________________________________________
694 Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Double_t yvertex,Double_t xvertex,Int_t pflag,AliITSclusterTable* table){
695 //function used to to find the clusters associated to the track
697 AliITSlayer &lay = fgLayers[layer];
698 Int_t * firstmod = new Int_t[fGeom->GetNlayers()];
699 for(Int_t i=0;i<fGeom->GetNlayers();i++){
700 firstmod[i]=fGeom->GetModuleIndex(i+1,1,1);
704 Float_t cx1,cx2,cy1,cy2;
705 FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3);
706 Int_t fun = FindIntersection(fCoef1,fCoef2,fCoef3,-(lay.GetR()*lay.GetR()),cx1,cy1,cx2,cy2);
708 Double_t fi1 =TMath::ATan2(cy1,cx1);
709 Double_t fi2 =TMath::ATan2(cy2,cx2);
710 fPhiEstimate = ChoosePoint(fi1,fi2,fPhic);
713 Double_t zed = TMath::Tan(fLambdac)*lay.GetR()+zvertex;
714 Double_t zed1 = TMath::Tan(fLambdac+lambdawindow)*lay.GetR()+zvertex;
715 Double_t zed2 = TMath::Tan(fLambdac-lambdawindow)*lay.GetR()+zvertex;
717 Double_t fi = fPhiEstimate;
718 Int_t nmod = lay.FindDetectorIndex(fi,zed)+firstmod[layer];
720 Int_t nm[8]={0,0,0,0,0,0,0,0};
721 nm[0] = lay.FindDetectorIndex(fi+phiwindow,zed)+firstmod[layer];
722 nm[1] = lay.FindDetectorIndex(fi-phiwindow,zed)+firstmod[layer];
723 nm[2] = lay.FindDetectorIndex(fi,zed1)+firstmod[layer];
724 nm[3] = lay.FindDetectorIndex(fi,zed2)+firstmod[layer];
725 nm[4] = lay.FindDetectorIndex(fi+phiwindow,zed1)+firstmod[layer];
726 nm[5] = lay.FindDetectorIndex(fi-phiwindow,zed1)+firstmod[layer];
727 nm[6] = lay.FindDetectorIndex(fi+phiwindow,zed2)+firstmod[layer];
728 nm[7] = lay.FindDetectorIndex(fi-phiwindow,zed2)+firstmod[layer];
733 TArrayI* array =(TArrayI*)table->GetListOfClusters(nmod);
734 TArrayI* list = new TArrayI(array->GetSize());
735 for(Int_t i=0;i<array->GetSize();i++){
736 Int_t in=(Int_t)array->At(i);
741 for(Int_t ii=0;ii<8;ii++){
742 if(nm[ii]!=nmod && nm[ii]>=0){
743 TArrayI* ar =(TArrayI*)table->GetListOfClusters(nm[ii]);
744 list->Set(list->GetSize()+ar->GetSize());
745 for(Int_t j=0;j<ar->GetSize();j++){
746 Int_t in=(Int_t)ar->At(j);
752 for(Int_t i=0;i<list->GetSize();i++){
753 Int_t index = (Int_t)list->At(i);
754 AliITSclusterV2* cllay = lay.GetCluster(index);
755 if(cllay==0) continue;
756 if(cllay->IsUsed()==1) continue;
759 Int_t modlay = cllay->GetDetectorIndex()+firstmod[layer];
760 GetClusterGCoordinates(cllay,modlay,x,y,z);
761 Double_t phi = TMath::ATan2(y,x);
762 Double_t lambda=TMath::ATan2(z-zvertex,TMath::Sqrt((x-xvertex)*(x-xvertex)+(y-yvertex)*(y-yvertex)));
765 if(TMath::Abs(fLambdac-lambda)<lambdawindow &&
766 TMath::Abs(fPhiEstimate-phi)<phiwindow){
769 if(trs->GetNumberOfClustersSA()==20){
773 trs->AddClusterSA(layer,index);
776 fPointc[0]=x;fPointc[1]=y;
786 //__________________________________________________________________
787 void AliITStrackerSA::GetClusterGCoordinates(AliITSclusterV2* cluster,Int_t module,Double_t& x, Double_t& y, Double_t& z){
789 //returns x,y,z of cluster in global coordinates
791 Double_t rot[9]; fGeom->GetRotMatrix(module,rot);
792 Int_t lay,lad,det; fGeom->GetModuleId(module,lay,lad,det);
793 Float_t tx,ty,tz; fGeom->GetTrans(lay,lad,det,tx,ty,tz);
795 Double_t alpha=TMath::ATan2(rot[1],rot[0])+TMath::Pi();
796 Double_t phi=TMath::Pi()/2+alpha;
797 if (lay==1) phi+=TMath::Pi();
799 Double_t cp=TMath::Cos(phi), sp=TMath::Sin(phi);
800 Double_t r=tx*cp+ty*sp;
802 x= r*cp - cluster->GetY()*sp;
803 y= r*sp + cluster->GetY()*cp;
809 //__________________________________________________________________
810 void AliITStrackerSA::GetClusterGErrors(AliITSclusterV2* cluster,Int_t module,Double_t& sigmax, Double_t& sigmay, Double_t& sigmaz){
812 //returns x,y,z of cluster in global coordinates
814 Double_t rot[9]; fGeom->GetRotMatrix(module,rot);
815 Int_t lay,lad,det; fGeom->GetModuleId(module,lay,lad,det);
817 Double_t alpha=TMath::ATan2(rot[1],rot[0])+TMath::Pi();
818 Double_t phi=TMath::Pi()/2+alpha;
819 if (lay==1) phi+=TMath::Pi();
821 Double_t cp=TMath::Cos(phi), sp=TMath::Sin(phi);
824 sigmax = TMath::Sqrt(sp*sp*cluster->GetSigmaY2());
825 sigmay = TMath::Sqrt(cp*cp*cluster->GetSigmaY2());
826 sigmaz = TMath::Sqrt(cluster->GetSigmaZ2());
829 //___________________________________________________________________
830 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){
832 //given (x,y) of three recpoints (in global coordinates)
833 //returns the parameters a,b,c of circonference x*x + y*y +a*x + b*y +c
835 Float_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
837 a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
838 b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
839 c = -x1*x1-y1*y1-a*x1-b*y1;
842 //__________________________________________________________________________
843 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){
845 //Finds the intersection between the circonference of the track and the circonference centered in (0,0) represented by one layer
846 //c2 is -rlayer*rlayer
850 Float_t aA = (b1*b1)/(a1*a1)+1;
851 Float_t bB = (-2*m*b1/(a1*a1));
852 Float_t cC = c2+(m*m)/(a1*a1);
853 if((bB*bB-4*aA*cC)<0) return 0;
855 y1 = (-bB+TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA);
856 y2 = (-bB-TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA);
857 x1 = (c2-c1-b1*y1)/a1;
858 x2 = (c2-c1-b1*y2)/a1;
862 //____________________________________________________________________
863 Double_t AliITStrackerSA::Curvature(Double_t x1,Double_t y1,Double_t
864 x2,Double_t y2,Double_t x3,Double_t y3){
866 //calculates the curvature of track
867 Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1);
869 Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den;
870 Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1);
871 Double_t c = -x1*x1-y1*y1-a*x1-b*y1;
874 if((a*a+b*b-4*c)<0) return 0;
875 Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.;
878 if((x1>0 && y1>0 && x1<xc)) rad*=-1;
879 if((x1<0 && y1>0 && x1<xc)) rad*=-1;
880 // if((x1<0 && y1<0 && x1<xc)) rad*=-1;
881 // if((x1>0 && y1<0 && x1<xc)) rad*=-1;
886 //____________________________________________________________________
887 Double_t AliITStrackerSA::ChoosePoint(Double_t p1, Double_t p2, Double_t pp){
889 //Returns the point closest to pp
891 Double_t diff1 = p1-pp;
892 Double_t diff2 = p2-pp;
894 if(TMath::Abs(diff1)<TMath::Abs(diff2)) fPhiEstimate=p1;
895 else fPhiEstimate=p2;
901 //_________________________________________________________________
902 AliITStrackV2* AliITStrackerSA::FindTrackLowChiSquare(TObjArray* tracklist, Int_t dim) const {
903 // returns track with lowes chi square
905 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(0);
909 Double_t * chi2 = new Double_t[dim];
910 Int_t * index = new Int_t[dim];
911 for(Int_t i=0;i<dim;i++){
912 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(i);
913 chi2[i]=trk->GetChi2();
917 Int_t w=0;Double_t value;
920 for(Int_t j=w+1;j<dim;j++){
933 AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(index[dim-1]);
940 //__________________________________________________________
941 Int_t AliITStrackerSA::FindLabel(Int_t l1, Int_t l2, Int_t l3, Int_t l4, Int_t l5, Int_t l6){
943 //function used to determine the track label
945 Int_t lb[6] = {l1,l2,l3,l4,l5,l6};
946 Int_t aa[6]={1,1,1,1,1,1};
949 Int_t k=0;Int_t w=0;Int_t num=6;
954 for(Int_t i=k+1;i<num;i++){
956 if(lb[k]==lb[i] && aa[k]!=0){
967 for(Int_t j=0;j<6;j++){
979 if(num==6) return lb[5];
983 //_____________________________________________________________________________
984 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,
985 Int_t gl12, Int_t gl13, Int_t gl14,Int_t gl15, Int_t gl16, Int_t gl17, Int_t gl18, Int_t numberofpoints){
988 //function used to assign label to the found track. If track is fake, the label is negative
990 Int_t lb0[6] = {gl1,gl2,gl3,gl4,gl5,gl6};
991 Int_t lb1[6] = {gl7,gl8,gl9,gl10,gl11,gl12};
992 Int_t lb2[6] = {gl13,gl14,gl15,gl16,gl17,gl18};
993 Int_t ll=FindLabel(lb0[0],lb0[1],lb0[2],lb0[3],lb0[4],lb0[5]);
994 Int_t lflag=0;Int_t num=6;
995 if(lb0[5]==-1 && lb1[5]==-1 && lb2[5]==-1) num=5;
997 for(Int_t i=0;i<num;i++){
998 if(lb0[i]==ll || lb1[i]==ll || lb2[i]==ll) lflag+=1;
1001 if(lflag>=numberofpoints) return ll;
1007 //_____________________________________________________________________________
1008 void AliITStrackerSA::SetWindowSizes(Int_t n, Double_t *phi, Double_t *lam){
1009 // Set sizes of the phi and lambda windows used for track finding
1011 if(phi){ // user defined values
1012 fPhiWin = new Double_t[fNloop];
1013 fLambdaWin = new Double_t[fNloop];
1014 for(Int_t k=0;k<fNloop;k++){
1016 fLambdaWin[k]=lam[k];
1019 else { // default values
1020 Double_t phid[46] = {0.001,0.0015,0.002,0.0023,0.0025,0.0027,0.003,
1021 0.0033,0.0035,0.0037,0.004,0.0043,0.0045,0.0047,
1022 0.005,0.0053,0.0055,
1023 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1024 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1025 0.01,0.0105,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1026 0.014,0.0145,0.015,0.0155,0.016,0.0165,0.017};
1027 Double_t lambdad[46] = {0.001,0.0015,0.002,0.0023,0.0025,0.0027,0.003,
1028 0.0033,0.0035,0.0037,0.004,0.0043,0.0045,0.0047,
1029 0.005,0.0053,0.0055,
1030 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077,
1031 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097,
1032 0.01,0.015,0.011,0.0115,0.012,0.0125,0.013,0.0135,
1033 0.014,0.0145,0.015,0.0155,0.016,0.016,0.016};
1036 Warning("SetWindowSizes","Number of loop forced to the default value %d",fNloop);
1039 fPhiWin = new Double_t[fNloop];
1040 fLambdaWin = new Double_t[fNloop];
1042 for(Int_t k=0;k<fNloop;k++){
1044 fLambdaWin[k]=lambdad[k];