1 /**************************************************************************
2 * Copyright(c) 1998-1999, 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 **************************************************************************/
17 //-------------------------------------------------------
18 // Implementation of the TPC tracker
20 // Origin: Marian Ivanov Marian.Ivanov@cern.ch
22 // AliTPC parallel tracker -
24 // run AliTPCFindClusters.C macro - clusters neccessary for tracker are founded
25 // run AliTPCFindTracksMI.C macro - to find tracks
26 // tracks are written to AliTPCtracks.root file
27 // for comparison also seeds are written to the same file - to special branch
28 //-------------------------------------------------------
33 #include "Riostream.h"
34 #include <TClonesArray.h>
36 #include <TObjArray.h>
39 #include "AliComplexCluster.h"
41 #include "AliESDkink.h"
43 #include "AliRunLoader.h"
44 #include "AliTPCClustersRow.h"
45 #include "AliTPCParam.h"
46 #include "AliTPCReconstructor.h"
47 #include "AliTPCclusterMI.h"
48 #include "AliTPCpolyTrack.h"
49 #include "AliTPCreco.h"
50 #include "AliTPCseed.h"
51 #include "AliTPCtrackerMI.h"
52 #include "TStopwatch.h"
53 #include "AliTPCReconstructor.h"
54 #include "AliESDkink.h"
56 #include "TTreeStream.h"
57 #include "AliAlignObj.h"
58 #include "AliTrackPointArray.h"
62 ClassImp(AliTPCtrackerMI)
65 class AliTPCFastMath {
68 static Double_t FastAsin(Double_t x);
70 static Double_t fgFastAsin[20000]; //lookup table for fast asin computation
73 Double_t AliTPCFastMath::fgFastAsin[20000];
74 AliTPCFastMath gAliTPCFastMath; // needed to fill the LUT
76 AliTPCFastMath::AliTPCFastMath(){
78 // initialized lookup table;
79 for (Int_t i=0;i<10000;i++){
80 fgFastAsin[2*i] = TMath::ASin(i/10000.);
81 fgFastAsin[2*i+1] = (TMath::ASin((i+1)/10000.)-fgFastAsin[2*i]);
85 Double_t AliTPCFastMath::FastAsin(Double_t x){
87 // return asin using lookup table
89 Int_t index = int(x*10000);
90 return fgFastAsin[2*index]+(x*10000.-index)*fgFastAsin[2*index+1];
93 Int_t index = int(x*10000);
94 return -(fgFastAsin[2*index]+(x*10000.-index)*fgFastAsin[2*index+1]);
100 Int_t AliTPCtrackerMI::UpdateTrack(AliTPCseed * track, Int_t accept){
102 //update track information using current cluster - track->fCurrentCluster
105 AliTPCclusterMI* c =track->fCurrentCluster;
106 if (accept>0) track->fCurrentClusterIndex1 |=0x8000; //sign not accepted clusters
108 UInt_t i = track->fCurrentClusterIndex1;
110 Int_t sec=(i&0xff000000)>>24;
111 //Int_t row = (i&0x00ff0000)>>16;
112 track->fRow=(i&0x00ff0000)>>16;
113 track->fSector = sec;
114 // Int_t index = i&0xFFFF;
115 if (sec>=fParam->GetNInnerSector()) track->fRow += fParam->GetNRowLow();
116 track->SetClusterIndex2(track->fRow, i);
117 //track->fFirstPoint = row;
118 //if ( track->fLastPoint<row) track->fLastPoint =row;
119 // if (track->fRow<0 || track->fRow>160) {
120 // printf("problem\n");
122 if (track->fFirstPoint>track->fRow)
123 track->fFirstPoint = track->fRow;
124 if (track->fLastPoint<track->fRow)
125 track->fLastPoint = track->fRow;
128 track->fClusterPointer[track->fRow] = c;
131 Float_t angle2 = track->GetSnp()*track->GetSnp();
132 angle2 = TMath::Sqrt(angle2/(1-angle2));
134 //SET NEW Track Point
138 AliTPCTrackerPoint &point =*(track->GetTrackPoint(track->fRow));
140 point.SetSigmaY(c->GetSigmaY2()/track->fCurrentSigmaY2);
141 point.SetSigmaZ(c->GetSigmaZ2()/track->fCurrentSigmaZ2);
142 point.SetErrY(sqrt(track->fErrorY2));
143 point.SetErrZ(sqrt(track->fErrorZ2));
145 point.SetX(track->GetX());
146 point.SetY(track->GetY());
147 point.SetZ(track->GetZ());
148 point.SetAngleY(angle2);
149 point.SetAngleZ(track->GetTgl());
150 if (point.fIsShared){
151 track->fErrorY2 *= 4;
152 track->fErrorZ2 *= 4;
156 Double_t chi2 = track->GetPredictedChi2(track->fCurrentCluster);
158 track->fErrorY2 *= 1.3;
159 track->fErrorY2 += 0.01;
160 track->fErrorZ2 *= 1.3;
161 track->fErrorZ2 += 0.005;
163 if (accept>0) return 0;
164 if (track->GetNumberOfClusters()%20==0){
165 // if (track->fHelixIn){
166 // TClonesArray & larr = *(track->fHelixIn);
167 // Int_t ihelix = larr.GetEntriesFast();
168 // new(larr[ihelix]) AliHelix(*track) ;
171 track->fNoCluster =0;
172 return track->Update(c,chi2,i);
177 Int_t AliTPCtrackerMI::AcceptCluster(AliTPCseed * seed, AliTPCclusterMI * cluster, Float_t factor,
178 Float_t cory, Float_t corz)
181 // decide according desired precision to accept given
182 // cluster for tracking
183 Double_t sy2=ErrY2(seed,cluster)*cory;
184 Double_t sz2=ErrZ2(seed,cluster)*corz;
185 //sy2=ErrY2(seed,cluster)*cory;
186 //sz2=ErrZ2(seed,cluster)*cory;
188 Double_t sdistancey2 = sy2+seed->GetSigmaY2();
189 Double_t sdistancez2 = sz2+seed->GetSigmaZ2();
191 Double_t rdistancey2 = (seed->fCurrentCluster->GetY()-seed->GetY())*
192 (seed->fCurrentCluster->GetY()-seed->GetY())/sdistancey2;
193 Double_t rdistancez2 = (seed->fCurrentCluster->GetZ()-seed->GetZ())*
194 (seed->fCurrentCluster->GetZ()-seed->GetZ())/sdistancez2;
196 Double_t rdistance2 = rdistancey2+rdistancez2;
199 if (rdistance2>16) return 3;
202 if ((rdistancey2>9.*factor || rdistancez2>9.*factor) && cluster->GetType()==0)
203 return 2; //suspisiouce - will be changed
205 if ((rdistancey2>6.25*factor || rdistancez2>6.25*factor) && cluster->GetType()>0)
206 // strict cut on overlaped cluster
207 return 2; //suspisiouce - will be changed
209 if ( (rdistancey2>1.*factor || rdistancez2>6.25*factor )
210 && cluster->GetType()<0){
220 //_____________________________________________________________________________
221 AliTPCtrackerMI::AliTPCtrackerMI(const AliTPCParam *par):
222 AliTracker(), fkNIS(par->GetNInnerSector()/2), fkNOS(par->GetNOuterSector()/2)
224 //---------------------------------------------------------------------
225 // The main TPC tracker constructor
226 //---------------------------------------------------------------------
227 fInnerSec=new AliTPCSector[fkNIS];
228 fOuterSec=new AliTPCSector[fkNOS];
231 for (i=0; i<fkNIS; i++) fInnerSec[i].Setup(par,0);
232 for (i=0; i<fkNOS; i++) fOuterSec[i].Setup(par,1);
239 Int_t nrowlow = par->GetNRowLow();
240 Int_t nrowup = par->GetNRowUp();
243 for (Int_t i=0;i<nrowlow;i++){
244 fXRow[i] = par->GetPadRowRadiiLow(i);
245 fPadLength[i]= par->GetPadPitchLength(0,i);
246 fYMax[i] = fXRow[i]*TMath::Tan(0.5*par->GetInnerAngle());
250 for (Int_t i=0;i<nrowup;i++){
251 fXRow[i+nrowlow] = par->GetPadRowRadiiUp(i);
252 fPadLength[i+nrowlow] = par->GetPadPitchLength(60,i);
253 fYMax[i+nrowlow] = fXRow[i+nrowlow]*TMath::Tan(0.5*par->GetOuterAngle());
264 fDebugStreamer = new TTreeSRedirector("TPCdebug.root");
266 //________________________________________________________________________
267 AliTPCtrackerMI::AliTPCtrackerMI(const AliTPCtrackerMI &t):
272 //------------------------------------
273 // dummy copy constructor
274 //------------------------------------------------------------------
276 AliTPCtrackerMI & AliTPCtrackerMI::operator=(const AliTPCtrackerMI& /*r*/){
277 //------------------------------
279 //--------------------------------------------------------------
282 //_____________________________________________________________________________
283 AliTPCtrackerMI::~AliTPCtrackerMI() {
284 //------------------------------------------------------------------
285 // TPC tracker destructor
286 //------------------------------------------------------------------
293 if (fDebugStreamer) delete fDebugStreamer;
296 void AliTPCtrackerMI::SetIO()
300 fInput = AliRunLoader::GetTreeR("TPC", kFALSE,AliConfig::GetDefaultEventFolderName());
302 fOutput = AliRunLoader::GetTreeT("TPC", kTRUE,AliConfig::GetDefaultEventFolderName());
304 AliTPCtrack *iotrack= new AliTPCtrack;
305 fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
311 void AliTPCtrackerMI::SetIO(TTree * input, TTree * output, AliESD * event)
327 AliTPCtrack *iotrack= new AliTPCtrack;
328 // iotrack->fHelixIn = new TClonesArray("AliHelix");
329 //iotrack->fHelixOut = new TClonesArray("AliHelix");
330 fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
333 if (output && (fDebug&2)){
334 //write the full seed information if specified in debug mode
336 fSeedTree = new TTree("Seeds","Seeds");
337 AliTPCseed * vseed = new AliTPCseed;
339 TClonesArray * arrtr = new TClonesArray("AliTPCTrackPoint",160);
340 arrtr->ExpandCreateFast(160);
341 TClonesArray * arre = new TClonesArray("AliTPCExactPoint",160);
343 vseed->fPoints = arrtr;
344 vseed->fEPoints = arre;
345 // vseed->fClusterPoints = arrcl;
346 fSeedTree->Branch("seeds","AliTPCseed",&vseed,32000,99);
349 fTreeDebug = new TTree("trackDebug","trackDebug");
350 TClonesArray * arrd = new TClonesArray("AliTPCTrackPoint2",0);
351 fTreeDebug->Branch("debug",&arrd,32000,99);
359 void AliTPCtrackerMI::FillESD(TObjArray* arr)
363 //fill esds using updated tracks
365 // write tracks to the event
366 // store index of the track
367 Int_t nseed=arr->GetEntriesFast();
368 //FindKinks(arr,fEvent);
369 for (Int_t i=0; i<nseed; i++) {
370 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
373 // pt->PropagateTo(fParam->GetInnerRadiusLow());
374 if (pt->GetKinkIndex(0)<=0){ //don't propagate daughter tracks
375 pt->PropagateTo(fParam->GetInnerRadiusLow());
378 if (( pt->GetPoints()[2]- pt->GetPoints()[0])>5 && pt->GetPoints()[3]>0.8){
380 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
381 iotrack.SetTPCPoints(pt->GetPoints());
382 iotrack.SetKinkIndexes(pt->GetKinkIndexes());
383 iotrack.SetV0Indexes(pt->GetV0Indexes());
384 // iotrack.SetTPCpid(pt->fTPCr);
385 //iotrack.SetTPCindex(i);
386 fEvent->AddTrack(&iotrack);
390 if ( (pt->GetNumberOfClusters()>70)&& (Float_t(pt->GetNumberOfClusters())/Float_t(pt->fNFoundable))>0.55) {
392 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
393 iotrack.SetTPCPoints(pt->GetPoints());
394 //iotrack.SetTPCindex(i);
395 iotrack.SetKinkIndexes(pt->GetKinkIndexes());
396 iotrack.SetV0Indexes(pt->GetV0Indexes());
397 // iotrack.SetTPCpid(pt->fTPCr);
398 fEvent->AddTrack(&iotrack);
402 // short tracks - maybe decays
404 if ( (pt->GetNumberOfClusters()>30) && (Float_t(pt->GetNumberOfClusters())/Float_t(pt->fNFoundable))>0.70) {
405 Int_t found,foundable,shared;
406 pt->GetClusterStatistic(0,60,found, foundable,shared,kFALSE);
407 if ( (found>20) && (pt->fNShared/float(pt->GetNumberOfClusters())<0.2)){
409 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
410 //iotrack.SetTPCindex(i);
411 iotrack.SetTPCPoints(pt->GetPoints());
412 iotrack.SetKinkIndexes(pt->GetKinkIndexes());
413 iotrack.SetV0Indexes(pt->GetV0Indexes());
414 //iotrack.SetTPCpid(pt->fTPCr);
415 fEvent->AddTrack(&iotrack);
420 if ( (pt->GetNumberOfClusters()>20) && (Float_t(pt->GetNumberOfClusters())/Float_t(pt->fNFoundable))>0.8) {
421 Int_t found,foundable,shared;
422 pt->GetClusterStatistic(0,60,found, foundable,shared,kFALSE);
423 if (found<20) continue;
424 if (pt->fNShared/float(pt->GetNumberOfClusters())>0.2) continue;
427 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
428 iotrack.SetTPCPoints(pt->GetPoints());
429 iotrack.SetKinkIndexes(pt->GetKinkIndexes());
430 iotrack.SetV0Indexes(pt->GetV0Indexes());
431 //iotrack.SetTPCpid(pt->fTPCr);
432 //iotrack.SetTPCindex(i);
433 fEvent->AddTrack(&iotrack);
436 // short tracks - secondaties
438 if ( (pt->GetNumberOfClusters()>30) ) {
439 Int_t found,foundable,shared;
440 pt->GetClusterStatistic(128,158,found, foundable,shared,kFALSE);
441 if ( (found>20) && (pt->fNShared/float(pt->GetNumberOfClusters())<0.2) &&float(found)/float(foundable)>0.8){
443 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
444 iotrack.SetTPCPoints(pt->GetPoints());
445 iotrack.SetKinkIndexes(pt->GetKinkIndexes());
446 iotrack.SetV0Indexes(pt->GetV0Indexes());
447 //iotrack.SetTPCpid(pt->fTPCr);
448 //iotrack.SetTPCindex(i);
449 fEvent->AddTrack(&iotrack);
454 if ( (pt->GetNumberOfClusters()>15)) {
455 Int_t found,foundable,shared;
456 pt->GetClusterStatistic(138,158,found, foundable,shared,kFALSE);
457 if (found<15) continue;
458 if (pt->fNShared/float(pt->GetNumberOfClusters())>0.2) continue;
459 if (float(found)/float(foundable)<0.8) continue;
462 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
463 iotrack.SetTPCPoints(pt->GetPoints());
464 iotrack.SetKinkIndexes(pt->GetKinkIndexes());
465 iotrack.SetV0Indexes(pt->GetV0Indexes());
466 // iotrack.SetTPCpid(pt->fTPCr);
467 //iotrack.SetTPCindex(i);
468 fEvent->AddTrack(&iotrack);
473 printf("Number of filled ESDs-\t%d\n",fEvent->GetNumberOfTracks());
476 void AliTPCtrackerMI::WriteTracks(TTree * tree)
479 // write tracks from seed array to selected tree
483 AliTPCtrack *iotrack= new AliTPCtrack;
484 fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
489 void AliTPCtrackerMI::WriteTracks()
492 // write tracks to the given output tree -
493 // output specified with SetIO routine
500 AliTPCtrack *iotrack= 0;
501 Int_t nseed=fSeeds->GetEntriesFast();
502 //for (Int_t i=0; i<nseed; i++) {
503 // iotrack= (AliTPCtrack*)fSeeds->UncheckedAt(i);
504 // if (iotrack) break;
506 //TBranch * br = fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
507 TBranch * br = fOutput->GetBranch("tracks");
508 br->SetAddress(&iotrack);
510 for (Int_t i=0; i<nseed; i++) {
511 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i);
513 AliTPCtrack * track = new AliTPCtrack(*pt);
516 // br->SetAddress(&iotrack);
521 //fOutput->GetDirectory()->cd();
527 //write the full seed information if specified in debug mode
529 AliTPCseed * vseed = new AliTPCseed;
531 TClonesArray * arrtr = new TClonesArray("AliTPCTrackPoint",160);
532 arrtr->ExpandCreateFast(160);
533 //TClonesArray * arrcl = new TClonesArray("AliTPCclusterMI",160);
534 //arrcl->ExpandCreateFast(160);
535 TClonesArray * arre = new TClonesArray("AliTPCExactPoint",160);
537 vseed->fPoints = arrtr;
538 vseed->fEPoints = arre;
539 // vseed->fClusterPoints = arrcl;
540 //TBranch * brseed = seedtree->Branch("seeds","AliTPCseed",&vseed,32000,99);
541 TBranch * brseed = fSeedTree->GetBranch("seeds");
543 Int_t nseed=fSeeds->GetEntriesFast();
545 for (Int_t i=0; i<nseed; i++) {
546 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i);
549 // pt->fClusterPoints = arrcl;
553 brseed->SetAddress(&vseed);
557 // pt->fClusterPoints = 0;
560 if (fTreeDebug) fTreeDebug->Write();
568 Double_t AliTPCtrackerMI::ErrY2(AliTPCseed* seed, AliTPCclusterMI * cl){
571 //seed->SetErrorY2(0.1);
573 //calculate look-up table at the beginning
574 static Bool_t ginit = kFALSE;
575 static Float_t gnoise1,gnoise2,gnoise3;
576 static Float_t ggg1[10000];
577 static Float_t ggg2[10000];
578 static Float_t ggg3[10000];
579 static Float_t glandau1[10000];
580 static Float_t glandau2[10000];
581 static Float_t glandau3[10000];
583 static Float_t gcor01[500];
584 static Float_t gcor02[500];
585 static Float_t gcorp[500];
590 for (Int_t i=1;i<500;i++){
591 Float_t rsigma = float(i)/100.;
592 gcor02[i] = TMath::Max(0.78 +TMath::Exp(7.4*(rsigma-1.2)),0.6);
593 gcor01[i] = TMath::Max(0.72 +TMath::Exp(3.36*(rsigma-1.2)),0.6);
594 gcorp[i] = TMath::Max(TMath::Power((rsigma+0.5),1.5),1.2);
598 for (Int_t i=3;i<10000;i++){
602 Float_t amp = float(i);
603 Float_t padlength =0.75;
604 gnoise1 = 0.0004/padlength;
605 Float_t nel = 0.268*amp;
606 Float_t nprim = 0.155*amp;
607 ggg1[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.001*nel/(padlength*padlength))/nel;
608 glandau1[i] = (2.+0.12*nprim)*0.5* (2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
609 if (glandau1[i]>1) glandau1[i]=1;
610 glandau1[i]*=padlength*padlength/12.;
614 gnoise2 = 0.0004/padlength;
617 ggg2[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
618 glandau2[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
619 if (glandau2[i]>1) glandau2[i]=1;
620 glandau2[i]*=padlength*padlength/12.;
625 gnoise3 = 0.0004/padlength;
628 ggg3[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
629 glandau3[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
630 if (glandau3[i]>1) glandau3[i]=1;
631 glandau3[i]*=padlength*padlength/12.;
639 Int_t amp = int(TMath::Abs(cl->GetQ()));
641 seed->SetErrorY2(1.);
645 Float_t z = TMath::Abs(fParam->GetZLength()-TMath::Abs(seed->GetZ()));
646 Int_t ctype = cl->GetType();
647 Float_t padlength= GetPadPitchLength(seed->fRow);
648 Float_t angle2 = seed->GetSnp()*seed->GetSnp();
649 angle2 = angle2/(1-angle2);
652 Int_t rsigmay = int(100.*cl->GetSigmaY2()/(seed->fCurrentSigmaY2));
655 if (fSectors==fInnerSec){
657 res = ggg1[amp]*z+glandau1[amp]*angle2;
658 if (ctype==0) res *= gcor01[rsigmay];
661 res*= gcorp[rsigmay];
667 res = ggg2[amp]*z+glandau2[amp]*angle2;
668 if (ctype==0) res *= gcor02[rsigmay];
671 res*= gcorp[rsigmay];
676 res = ggg3[amp]*z+glandau3[amp]*angle2;
677 if (ctype==0) res *= gcor02[rsigmay];
680 res*= gcorp[rsigmay];
687 res*=2.4; // overestimate error 2 times
694 seed->SetErrorY2(res);
702 Double_t AliTPCtrackerMI::ErrZ2(AliTPCseed* seed, AliTPCclusterMI * cl){
705 //seed->SetErrorY2(0.1);
707 //calculate look-up table at the beginning
708 static Bool_t ginit = kFALSE;
709 static Float_t gnoise1,gnoise2,gnoise3;
710 static Float_t ggg1[10000];
711 static Float_t ggg2[10000];
712 static Float_t ggg3[10000];
713 static Float_t glandau1[10000];
714 static Float_t glandau2[10000];
715 static Float_t glandau3[10000];
717 static Float_t gcor01[1000];
718 static Float_t gcor02[1000];
719 static Float_t gcorp[1000];
724 for (Int_t i=1;i<1000;i++){
725 Float_t rsigma = float(i)/100.;
726 gcor02[i] = TMath::Max(0.81 +TMath::Exp(6.8*(rsigma-1.2)),0.6);
727 gcor01[i] = TMath::Max(0.72 +TMath::Exp(2.04*(rsigma-1.2)),0.6);
728 gcorp[i] = TMath::Max(TMath::Power((rsigma+0.5),1.5),1.2);
732 for (Int_t i=3;i<10000;i++){
736 Float_t amp = float(i);
737 Float_t padlength =0.75;
738 gnoise1 = 0.0004/padlength;
739 Float_t nel = 0.268*amp;
740 Float_t nprim = 0.155*amp;
741 ggg1[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.001*nel/(padlength*padlength))/nel;
742 glandau1[i] = (2.+0.12*nprim)*0.5* (2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
743 if (glandau1[i]>1) glandau1[i]=1;
744 glandau1[i]*=padlength*padlength/12.;
748 gnoise2 = 0.0004/padlength;
751 ggg2[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
752 glandau2[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
753 if (glandau2[i]>1) glandau2[i]=1;
754 glandau2[i]*=padlength*padlength/12.;
759 gnoise3 = 0.0004/padlength;
762 ggg3[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
763 glandau3[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
764 if (glandau3[i]>1) glandau3[i]=1;
765 glandau3[i]*=padlength*padlength/12.;
773 Int_t amp = int(TMath::Abs(cl->GetQ()));
775 seed->SetErrorY2(1.);
779 Float_t z = TMath::Abs(fParam->GetZLength()-TMath::Abs(seed->GetZ()));
780 Int_t ctype = cl->GetType();
781 Float_t padlength= GetPadPitchLength(seed->fRow);
783 Float_t angle2 = seed->GetSnp()*seed->GetSnp();
784 // if (angle2<0.6) angle2 = 0.6;
785 angle2 = seed->GetTgl()*seed->GetTgl()*(1+angle2/(1-angle2));
788 Int_t rsigmaz = int(100.*cl->GetSigmaZ2()/(seed->fCurrentSigmaZ2));
791 if (fSectors==fInnerSec){
793 res = ggg1[amp]*z+glandau1[amp]*angle2;
794 if (ctype==0) res *= gcor01[rsigmaz];
797 res*= gcorp[rsigmaz];
803 res = ggg2[amp]*z+glandau2[amp]*angle2;
804 if (ctype==0) res *= gcor02[rsigmaz];
807 res*= gcorp[rsigmaz];
812 res = ggg3[amp]*z+glandau3[amp]*angle2;
813 if (ctype==0) res *= gcor02[rsigmaz];
816 res*= gcorp[rsigmaz];
825 if ((ctype<0) &&<70){
833 seed->SetErrorZ2(res);
840 Double_t AliTPCtrackerMI::ErrZ2(AliTPCseed* seed, AliTPCclusterMI * cl){
843 //seed->SetErrorZ2(0.1);
847 Float_t z = TMath::Abs(fParam->GetZLength()-TMath::Abs(seed->GetZ()));
849 Float_t rsigmaz = cl->GetSigmaZ2()/(seed->fCurrentSigmaZ2);
850 Int_t ctype = cl->GetType();
851 Float_t amp = TMath::Abs(cl->GetQ());
856 Float_t landau=2 ; //landau fluctuation part
857 Float_t gg=2; // gg fluctuation part
858 Float_t padlength= GetPadPitchLength(seed->GetX());
860 if (fSectors==fInnerSec){
861 snoise2 = 0.0004/padlength;
864 gg = (2+0.001*nel/(padlength*padlength))/nel;
865 landau = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
866 if (landau>1) landau=1;
869 snoise2 = 0.0004/padlength;
872 gg = (2+0.0008*nel/(padlength*padlength))/nel;
873 landau = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
874 if (landau>1) landau=1;
876 Float_t sdiff = gg*fParam->GetDiffT()*fParam->GetDiffT()*z;
879 Float_t angle2 = seed->GetSnp()*seed->GetSnp();
880 angle2 = TMath::Sqrt((1-angle2));
881 if (angle2<0.6) angle2 = 0.6;
884 Float_t angle = seed->GetTgl()/angle2;
885 Float_t angular = landau*angle*angle*padlength*padlength/12.;
886 Float_t res = sdiff + angular;
889 if ((ctype==0) && (fSectors ==fOuterSec))
890 res *= 0.81 +TMath::Exp(6.8*(rsigmaz-1.2));
892 if ((ctype==0) && (fSectors ==fInnerSec))
893 res *= 0.72 +TMath::Exp(2.04*(rsigmaz-1.2));
897 res*= TMath::Power(rsigmaz+0.5,1.5); //0.31+0.147*ctype;
903 if ((ctype<0) &&<70){
911 seed->SetErrorZ2(res);
919 void AliTPCtrackerMI::RotateToLocal(AliTPCseed *seed)
921 //rotate to track "local coordinata
922 Float_t x = seed->GetX();
923 Float_t y = seed->GetY();
924 Float_t ymax = x*TMath::Tan(0.5*fSectors->GetAlpha());
927 seed->fRelativeSector= (seed->fRelativeSector+1) % fN;
928 if (!seed->Rotate(fSectors->GetAlpha()))
930 } else if (y <-ymax) {
931 seed->fRelativeSector= (seed->fRelativeSector-1+fN) % fN;
932 if (!seed->Rotate(-fSectors->GetAlpha()))
940 //_____________________________________________________________________________
941 Double_t AliTPCtrackerMI::F1old(Double_t x1,Double_t y1,
942 Double_t x2,Double_t y2,
943 Double_t x3,Double_t y3)
945 //-----------------------------------------------------------------
946 // Initial approximation of the track curvature
947 //-----------------------------------------------------------------
948 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
949 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
950 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
951 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
952 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
954 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
955 if ( xr*xr+yr*yr<=0.00000000000001) return 100;
956 return -xr*yr/sqrt(xr*xr+yr*yr);
961 //_____________________________________________________________________________
962 Double_t AliTPCtrackerMI::F1(Double_t x1,Double_t y1,
963 Double_t x2,Double_t y2,
964 Double_t x3,Double_t y3)
966 //-----------------------------------------------------------------
967 // Initial approximation of the track curvature
968 //-----------------------------------------------------------------
974 Double_t det = x3*y2-x2*y3;
979 Double_t u = 0.5* (x2*(x2-x3)+y2*(y2-y3))/det;
980 Double_t x0 = x3*0.5-y3*u;
981 Double_t y0 = y3*0.5+x3*u;
982 Double_t c2 = 1/TMath::Sqrt(x0*x0+y0*y0);
988 Double_t AliTPCtrackerMI::F2(Double_t x1,Double_t y1,
989 Double_t x2,Double_t y2,
990 Double_t x3,Double_t y3)
992 //-----------------------------------------------------------------
993 // Initial approximation of the track curvature
994 //-----------------------------------------------------------------
1000 Double_t det = x3*y2-x2*y3;
1005 Double_t u = 0.5* (x2*(x2-x3)+y2*(y2-y3))/det;
1006 Double_t x0 = x3*0.5-y3*u;
1007 Double_t y0 = y3*0.5+x3*u;
1008 Double_t c2 = 1/TMath::Sqrt(x0*x0+y0*y0);
1017 //_____________________________________________________________________________
1018 Double_t AliTPCtrackerMI::F2old(Double_t x1,Double_t y1,
1019 Double_t x2,Double_t y2,
1020 Double_t x3,Double_t y3)
1022 //-----------------------------------------------------------------
1023 // Initial approximation of the track curvature times center of curvature
1024 //-----------------------------------------------------------------
1025 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
1026 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
1027 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
1028 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
1029 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
1031 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
1033 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
1036 //_____________________________________________________________________________
1037 Double_t AliTPCtrackerMI::F3(Double_t x1,Double_t y1,
1038 Double_t x2,Double_t y2,
1039 Double_t z1,Double_t z2)
1041 //-----------------------------------------------------------------
1042 // Initial approximation of the tangent of the track dip angle
1043 //-----------------------------------------------------------------
1044 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
1048 Double_t AliTPCtrackerMI::F3n(Double_t x1,Double_t y1,
1049 Double_t x2,Double_t y2,
1050 Double_t z1,Double_t z2, Double_t c)
1052 //-----------------------------------------------------------------
1053 // Initial approximation of the tangent of the track dip angle
1054 //-----------------------------------------------------------------
1058 //angle1 = (z1-z2)*c/(TMath::ASin(c*x1-ni)-TMath::ASin(c*x2-ni));
1060 Double_t d = TMath::Sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
1061 if (TMath::Abs(d*c*0.5)>1) return 0;
1062 // Double_t angle2 = TMath::ASin(d*c*0.5);
1063 // Double_t angle2 = AliTPCFastMath::FastAsin(d*c*0.5);
1064 Double_t angle2 = (d*c*0.5>0.1)? TMath::ASin(d*c*0.5): AliTPCFastMath::FastAsin(d*c*0.5);
1066 angle2 = (z1-z2)*c/(angle2*2.);
1070 Bool_t AliTPCtrackerMI::GetProlongation(Double_t x1, Double_t x2, Double_t x[5], Double_t &y, Double_t &z)
1071 {//-----------------------------------------------------------------
1072 // This function find proloncation of a track to a reference plane x=x2.
1073 //-----------------------------------------------------------------
1077 if (TMath::Abs(x[4]*x1 - x[2]) >= 0.999) {
1081 Double_t c1=x[4]*x1 - x[2], r1=sqrt(1.- c1*c1);
1082 Double_t c2=x[4]*x2 - x[2], r2=sqrt(1.- c2*c2);
1086 Double_t dy = dx*(c1+c2)/(r1+r2);
1089 Double_t delta = x[4]*dx*(c1+c2)/(c1*r2 + c2*r1);
1091 if (TMath::Abs(delta)>0.01){
1092 dz = x[3]*TMath::ASin(delta)/x[4];
1094 dz = x[3]*AliTPCFastMath::FastAsin(delta)/x[4];
1097 //dz = x[3]*AliTPCFastMath::FastAsin(delta)/x[4];
1105 Int_t AliTPCtrackerMI::LoadClusters (TTree *tree)
1110 return LoadClusters();
1113 Int_t AliTPCtrackerMI::LoadClusters()
1116 // load clusters to the memory
1117 AliTPCClustersRow *clrow= new AliTPCClustersRow;
1118 clrow->SetClass("AliTPCclusterMI");
1120 clrow->GetArray()->ExpandCreateFast(10000);
1122 // TTree * tree = fClustersArray.GetTree();
1124 TTree * tree = fInput;
1125 TBranch * br = tree->GetBranch("Segment");
1126 br->SetAddress(&clrow);
1128 Int_t j=Int_t(tree->GetEntries());
1129 for (Int_t i=0; i<j; i++) {
1133 fParam->AdjustSectorRow(clrow->GetID(),sec,row);
1135 AliTPCRow * tpcrow=0;
1138 tpcrow = &(fInnerSec[sec%fkNIS][row]);
1142 tpcrow = &(fOuterSec[(sec-fkNIS*2)%fkNOS][row]);
1143 left = (sec-fkNIS*2)/fkNOS;
1146 tpcrow->fN1 = clrow->GetArray()->GetEntriesFast();
1147 tpcrow->fClusters1 = new AliTPCclusterMI[tpcrow->fN1];
1148 for (Int_t i=0;i<tpcrow->fN1;i++)
1149 tpcrow->fClusters1[i] = *(AliTPCclusterMI*)(clrow->GetArray()->At(i));
1152 tpcrow->fN2 = clrow->GetArray()->GetEntriesFast();
1153 tpcrow->fClusters2 = new AliTPCclusterMI[tpcrow->fN2];
1154 for (Int_t i=0;i<tpcrow->fN2;i++)
1155 tpcrow->fClusters2[i] = *(AliTPCclusterMI*)(clrow->GetArray()->At(i));
1166 void AliTPCtrackerMI::UnloadClusters()
1169 // unload clusters from the memory
1171 Int_t nrows = fOuterSec->GetNRows();
1172 for (Int_t sec = 0;sec<fkNOS;sec++)
1173 for (Int_t row = 0;row<nrows;row++){
1174 AliTPCRow* tpcrow = &(fOuterSec[sec%fkNOS][row]);
1176 // if (tpcrow->fClusters1) delete []tpcrow->fClusters1;
1177 // if (tpcrow->fClusters2) delete []tpcrow->fClusters2;
1179 tpcrow->ResetClusters();
1182 nrows = fInnerSec->GetNRows();
1183 for (Int_t sec = 0;sec<fkNIS;sec++)
1184 for (Int_t row = 0;row<nrows;row++){
1185 AliTPCRow* tpcrow = &(fInnerSec[sec%fkNIS][row]);
1187 // if (tpcrow->fClusters1) delete []tpcrow->fClusters1;
1188 //if (tpcrow->fClusters2) delete []tpcrow->fClusters2;
1190 tpcrow->ResetClusters();
1197 //_____________________________________________________________________________
1198 Int_t AliTPCtrackerMI::LoadOuterSectors() {
1199 //-----------------------------------------------------------------
1200 // This function fills outer TPC sectors with clusters.
1201 //-----------------------------------------------------------------
1202 Int_t nrows = fOuterSec->GetNRows();
1204 for (Int_t sec = 0;sec<fkNOS;sec++)
1205 for (Int_t row = 0;row<nrows;row++){
1206 AliTPCRow* tpcrow = &(fOuterSec[sec%fkNOS][row]);
1207 Int_t sec2 = sec+2*fkNIS;
1209 Int_t ncl = tpcrow->fN1;
1211 AliTPCclusterMI *c= &(tpcrow->fClusters1[ncl]);
1212 index=(((sec2<<8)+row)<<16)+ncl;
1213 tpcrow->InsertCluster(c,index);
1218 AliTPCclusterMI *c= &(tpcrow->fClusters2[ncl]);
1219 index=((((sec2+fkNOS)<<8)+row)<<16)+ncl;
1220 tpcrow->InsertCluster(c,index);
1223 // write indexes for fast acces
1225 for (Int_t i=0;i<510;i++)
1226 tpcrow->fFastCluster[i]=-1;
1227 for (Int_t i=0;i<tpcrow->GetN();i++){
1228 Int_t zi = Int_t((*tpcrow)[i]->GetZ()+255.);
1229 tpcrow->fFastCluster[zi]=i; // write index
1232 for (Int_t i=0;i<510;i++){
1233 if (tpcrow->fFastCluster[i]<0)
1234 tpcrow->fFastCluster[i] = last;
1236 last = tpcrow->fFastCluster[i];
1245 //_____________________________________________________________________________
1246 Int_t AliTPCtrackerMI::LoadInnerSectors() {
1247 //-----------------------------------------------------------------
1248 // This function fills inner TPC sectors with clusters.
1249 //-----------------------------------------------------------------
1250 Int_t nrows = fInnerSec->GetNRows();
1252 for (Int_t sec = 0;sec<fkNIS;sec++)
1253 for (Int_t row = 0;row<nrows;row++){
1254 AliTPCRow* tpcrow = &(fInnerSec[sec%fkNIS][row]);
1257 Int_t ncl = tpcrow->fN1;
1259 AliTPCclusterMI *c= &(tpcrow->fClusters1[ncl]);
1260 index=(((sec<<8)+row)<<16)+ncl;
1261 tpcrow->InsertCluster(c,index);
1266 AliTPCclusterMI *c= &(tpcrow->fClusters2[ncl]);
1267 index=((((sec+fkNIS)<<8)+row)<<16)+ncl;
1268 tpcrow->InsertCluster(c,index);
1271 // write indexes for fast acces
1273 for (Int_t i=0;i<510;i++)
1274 tpcrow->fFastCluster[i]=-1;
1275 for (Int_t i=0;i<tpcrow->GetN();i++){
1276 Int_t zi = Int_t((*tpcrow)[i]->GetZ()+255.);
1277 tpcrow->fFastCluster[zi]=i; // write index
1280 for (Int_t i=0;i<510;i++){
1281 if (tpcrow->fFastCluster[i]<0)
1282 tpcrow->fFastCluster[i] = last;
1284 last = tpcrow->fFastCluster[i];
1296 //_________________________________________________________________________
1297 AliTPCclusterMI *AliTPCtrackerMI::GetClusterMI(Int_t index) const {
1298 //--------------------------------------------------------------------
1299 // Return pointer to a given cluster
1300 //--------------------------------------------------------------------
1301 Int_t sec=(index&0xff000000)>>24;
1302 Int_t row=(index&0x00ff0000)>>16;
1303 Int_t ncl=(index&0x00007fff)>>00;
1305 const AliTPCRow * tpcrow=0;
1306 AliTPCclusterMI * clrow =0;
1309 tpcrow = &(fInnerSec[sec%fkNIS][row]);
1310 if (tpcrow==0) return 0;
1313 if (tpcrow->fN1<=ncl) return 0;
1314 clrow = tpcrow->fClusters1;
1317 if (tpcrow->fN2<=ncl) return 0;
1318 clrow = tpcrow->fClusters2;
1322 tpcrow = &(fOuterSec[(sec-fkNIS*2)%fkNOS][row]);
1323 if (tpcrow==0) return 0;
1325 if (sec-2*fkNIS<fkNOS) {
1326 if (tpcrow->fN1<=ncl) return 0;
1327 clrow = tpcrow->fClusters1;
1330 if (tpcrow->fN2<=ncl) return 0;
1331 clrow = tpcrow->fClusters2;
1335 return &(clrow[ncl]);
1341 Int_t AliTPCtrackerMI::FollowToNext(AliTPCseed& t, Int_t nr) {
1342 //-----------------------------------------------------------------
1343 // This function tries to find a track prolongation to next pad row
1344 //-----------------------------------------------------------------
1346 Double_t x= GetXrow(nr), ymax=GetMaxY(nr);
1347 AliTPCclusterMI *cl=0;
1348 Int_t tpcindex= t.GetClusterIndex2(nr);
1350 // update current shape info every 5 pad-row
1351 // if ( (nr%5==0) || t.GetNumberOfClusters()<2 || (t.fCurrentSigmaY2<0.0001) ){
1355 if (fIteration>0 && tpcindex>=-1){ //if we have already clusters
1357 if (tpcindex==-1) return 0; //track in dead zone
1359 cl = t.fClusterPointer[nr];
1360 if ( (cl==0) ) cl = GetClusterMI(tpcindex);
1361 t.fCurrentClusterIndex1 = tpcindex;
1364 Int_t relativesector = ((tpcindex&0xff000000)>>24)%18; // if previously accepted cluster in different sector
1365 Float_t angle = relativesector*fSectors->GetAlpha()+fSectors->GetAlphaShift();
1367 if (angle<-TMath::Pi()) angle += 2*TMath::Pi();
1368 if (angle>=TMath::Pi()) angle -= 2*TMath::Pi();
1370 if (TMath::Abs(angle-t.GetAlpha())>0.001){
1371 Double_t rotation = angle-t.GetAlpha();
1372 t.fRelativeSector= relativesector;
1377 t.fCurrentCluster = cl;
1379 Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1380 if ((tpcindex&0x8000)==0) accept =0;
1382 //if founded cluster is acceptible
1383 if (cl->IsUsed(11)) { // id cluster is shared inrease uncertainty
1390 UpdateTrack(&t,accept);
1395 if (fIteration>1) {t.fNFoundable++; return 0;} // not look for new cluster during refitting
1398 if (TMath::Abs(t.GetSnp())>0.95 || TMath::Abs(x*t.GetC()-t.GetEta())>0.95) return 0;
1399 Double_t y=t.GetYat(x);
1400 if (TMath::Abs(y)>ymax){
1402 t.fRelativeSector= (t.fRelativeSector+1) % fN;
1403 if (!t.Rotate(fSectors->GetAlpha()))
1405 } else if (y <-ymax) {
1406 t.fRelativeSector= (t.fRelativeSector-1+fN) % fN;
1407 if (!t.Rotate(-fSectors->GetAlpha()))
1413 if (!t.PropagateTo(x)) {
1414 if (fIteration==0) t.fRemoval = 10;
1418 Double_t z=t.GetZ();
1420 const AliTPCRow &krow=GetRow(t.fRelativeSector,nr);
1421 if ( (t.GetSigmaY2()<0) || t.GetSigmaZ2()<0) return 0;
1423 Double_t roadz = 1.;
1425 if (TMath::Abs(TMath::Abs(y)-ymax)<krow.fDeadZone){
1427 t.SetClusterIndex2(nr,-1);
1432 if (TMath::Abs(z)<(AliTPCReconstructor::GetCtgRange()*x+10) && TMath::Abs(z)<fParam->GetZLength() ) t.fNFoundable++;
1438 // cl = krow.FindNearest2(y+10.,z,roady,roadz,index);
1439 cl = krow.FindNearest2(y,z,roady,roadz,index);
1440 if (cl) t.fCurrentClusterIndex1 = krow.GetIndex(index);
1443 t.fCurrentCluster = cl;
1445 if (fIteration==2&&cl->IsUsed(10)) return 0;
1446 Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1447 if (fIteration==2&&cl->IsUsed(11)) {
1454 if (t.fCurrentCluster->IsUsed(10)){
1459 if (t.fNShared>0.7*t.GetNumberOfClusters()) {
1465 if (accept<3) UpdateTrack(&t,accept);
1468 if ( fIteration==0 && t.fNFoundable*0.5 > t.GetNumberOfClusters()) t.fRemoval=10;
1474 Int_t AliTPCtrackerMI::FollowToNextFast(AliTPCseed& t, Int_t nr) {
1475 //-----------------------------------------------------------------
1476 // This function tries to find a track prolongation to next pad row
1477 //-----------------------------------------------------------------
1479 Double_t x= GetXrow(nr), ymax=GetMaxY(nr);
1481 if (!t.GetProlongation(x,y,z)) {
1487 if (TMath::Abs(y)>ymax){
1490 t.fRelativeSector= (t.fRelativeSector+1) % fN;
1491 if (!t.Rotate(fSectors->GetAlpha()))
1493 } else if (y <-ymax) {
1494 t.fRelativeSector= (t.fRelativeSector-1+fN) % fN;
1495 if (!t.Rotate(-fSectors->GetAlpha()))
1498 if (!t.PropagateTo(x)) {
1501 t.GetProlongation(x,y,z);
1504 // update current shape info every 3 pad-row
1505 if ( (nr%6==0) || t.GetNumberOfClusters()<2 || (t.fCurrentSigmaY2<0.0001) ){
1506 // t.fCurrentSigmaY = GetSigmaY(&t);
1507 //t.fCurrentSigmaZ = GetSigmaZ(&t);
1511 AliTPCclusterMI *cl=0;
1516 const AliTPCRow &krow=GetRow(t.fRelativeSector,nr);
1517 if ( (t.GetSigmaY2()<0) || t.GetSigmaZ2()<0) return 0;
1519 Double_t roadz = 1.;
1522 if (TMath::Abs(TMath::Abs(y)-ymax)<krow.fDeadZone){
1524 t.SetClusterIndex2(row,-1);
1529 if (TMath::Abs(z)>(AliTPCReconstructor::GetCtgRange()*x+10)) t.SetClusterIndex2(row,-1);
1533 if ((cl==0)&&(krow)) {
1534 // cl = krow.FindNearest2(y+10,z,roady,roadz,index);
1535 cl = krow.FindNearest2(y,z,roady,roadz,index);
1537 if (cl) t.fCurrentClusterIndex1 = krow.GetIndex(index);
1541 t.fCurrentCluster = cl;
1542 // Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1544 t.SetClusterIndex2(row,index);
1545 t.fClusterPointer[row] = cl;
1553 //_________________________________________________________________________
1554 Bool_t AliTPCtrackerMI::GetTrackPoint(Int_t index, AliTrackPoint &p ) const
1556 // Get track space point by index
1557 // return false in case the cluster doesn't exist
1558 AliTPCclusterMI *cl = GetClusterMI(index);
1559 if (!cl) return kFALSE;
1560 Int_t sector = (index&0xff000000)>>24;
1561 Int_t row = (index&0x00ff0000)>>16;
1563 xyz[0] = fParam->GetPadRowRadii(sector,row);
1564 xyz[1] = cl->GetY();
1565 xyz[2] = cl->GetZ();
1567 fParam->AdjustCosSin(sector,cos,sin);
1568 Float_t x = cos*xyz[0]-sin*xyz[1];
1569 Float_t y = cos*xyz[1]+sin*xyz[0];
1571 Float_t sigmaY2 = 0.027*cl->GetSigmaY2();
1572 if (sector < fParam->GetNInnerSector()) sigmaY2 *= 2.07;
1573 Float_t sigmaZ2 = 0.066*cl->GetSigmaZ2();
1574 if (sector < fParam->GetNInnerSector()) sigmaZ2 *= 1.77;
1575 cov[0] = sin*sin*sigmaY2;
1576 cov[1] = -sin*cos*sigmaY2;
1578 cov[3] = cos*cos*sigmaY2;
1581 p.SetXYZ(x,y,xyz[2],cov);
1582 AliAlignObj::ELayerID iLayer;
1584 if (sector < fParam->GetNInnerSector()) {
1585 iLayer = AliAlignObj::kTPC1;
1589 iLayer = AliAlignObj::kTPC2;
1590 idet = sector - fParam->GetNInnerSector();
1592 UShort_t volid = AliAlignObj::LayerToVolUID(iLayer,idet);
1593 p.SetVolumeID(volid);
1599 Int_t AliTPCtrackerMI::UpdateClusters(AliTPCseed& t, Int_t nr) {
1600 //-----------------------------------------------------------------
1601 // This function tries to find a track prolongation to next pad row
1602 //-----------------------------------------------------------------
1603 t.fCurrentCluster = 0;
1604 t.fCurrentClusterIndex1 = 0;
1606 Double_t xt=t.GetX();
1607 Int_t row = GetRowNumber(xt)-1;
1608 Double_t ymax= GetMaxY(nr);
1610 if (row < nr) return 1; // don't prolongate if not information until now -
1611 if (TMath::Abs(t.GetSnp())>0.9 && t.GetNumberOfClusters()>40. && fIteration!=2) {
1613 return 0; // not prolongate strongly inclined tracks
1615 if (TMath::Abs(t.GetSnp())>0.95) {
1617 return 0; // not prolongate strongly inclined tracks
1620 Double_t x= GetXrow(nr);
1622 //t.PropagateTo(x+0.02);
1623 //t.PropagateTo(x+0.01);
1624 if (!t.PropagateTo(x)){
1631 if (TMath::Abs(y)>ymax){
1633 t.fRelativeSector= (t.fRelativeSector+1) % fN;
1634 if (!t.Rotate(fSectors->GetAlpha()))
1636 } else if (y <-ymax) {
1637 t.fRelativeSector= (t.fRelativeSector-1+fN) % fN;
1638 if (!t.Rotate(-fSectors->GetAlpha()))
1641 // if (!t.PropagateTo(x)){
1649 AliTPCRow &krow=GetRow(t.fRelativeSector,nr);
1651 if (TMath::Abs(TMath::Abs(y)-ymax)<krow.fDeadZone){
1653 t.SetClusterIndex2(nr,-1);
1658 if (TMath::Abs(t.GetZ())<(AliTPCReconstructor::GetCtgRange()*t.GetX()+10)) t.fNFoundable++;
1664 if ( (nr%6==0) || t.GetNumberOfClusters()<2){
1665 // t.fCurrentSigmaY = GetSigmaY(&t);
1666 //t.fCurrentSigmaZ = GetSigmaZ(&t);
1670 AliTPCclusterMI *cl=0;
1673 Double_t roady = 1.;
1674 Double_t roadz = 1.;
1678 index = t.GetClusterIndex2(nr);
1679 if ( (index>0) && (index&0x8000)==0){
1680 cl = t.fClusterPointer[nr];
1681 if ( (cl==0) && (index>0)) cl = GetClusterMI(index);
1682 t.fCurrentClusterIndex1 = index;
1684 t.fCurrentCluster = cl;
1690 // if (index<0) return 0;
1691 UInt_t uindex = TMath::Abs(index);
1694 //cl = krow.FindNearest2(y+10,z,roady,roadz,uindex);
1695 cl = krow.FindNearest2(y,z,roady,roadz,uindex);
1698 if (cl) t.fCurrentClusterIndex1 = krow.GetIndex(uindex);
1699 t.fCurrentCluster = cl;
1705 Int_t AliTPCtrackerMI::FollowToNextCluster(AliTPCseed & t, Int_t nr) {
1706 //-----------------------------------------------------------------
1707 // This function tries to find a track prolongation to next pad row
1708 //-----------------------------------------------------------------
1710 //update error according neighborhoud
1712 if (t.fCurrentCluster) {
1714 Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1716 if (t.fCurrentCluster->IsUsed(10)){
1722 if (t.fNShared>0.7*t.GetNumberOfClusters()) {
1727 if (fIteration>0) accept = 0;
1728 if (accept<3) UpdateTrack(&t,accept);
1732 if ( ( (t.GetSigmaY2()+t.GetSigmaZ2())>0.16)&& t.GetNumberOfClusters()>18) t.fRemoval=10;
1733 if ( t.GetChi2()/t.GetNumberOfClusters()>6 &&t.GetNumberOfClusters()>18) t.fRemoval=10;
1735 if (( (t.fNFoundable*0.5 > t.GetNumberOfClusters()) || t.fNoCluster>15)) t.fRemoval=10;
1743 //_____________________________________________________________________________
1744 Int_t AliTPCtrackerMI::FollowProlongation(AliTPCseed& t, Int_t rf, Int_t step) {
1745 //-----------------------------------------------------------------
1746 // This function tries to find a track prolongation.
1747 //-----------------------------------------------------------------
1748 Double_t xt=t.GetX();
1750 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
1751 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1752 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1754 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
1756 Int_t first = GetRowNumber(xt)-1;
1757 for (Int_t nr= first; nr>=rf; nr-=step) {
1759 if (t.GetKinkIndexes()[0]>0){
1760 for (Int_t i=0;i<3;i++){
1761 Int_t index = t.GetKinkIndexes()[i];
1762 if (index==0) break;
1763 if (index<0) continue;
1765 AliESDkink * kink = fEvent->GetKink(index-1);
1767 printf("PROBLEM\n");
1770 Int_t kinkrow = kink->GetTPCRow0()+2+Int_t(0.5/(0.05+kink->GetAngle(2)));
1772 AliExternalTrackParam paramd(t);
1773 kink->SetDaughter(paramd);
1774 kink->SetStatus(2,5);
1781 if (nr==80) t.UpdateReference();
1782 if (nr<fInnerSec->GetNRows())
1783 fSectors = fInnerSec;
1785 fSectors = fOuterSec;
1786 if (FollowToNext(t,nr)==0)
1795 //_____________________________________________________________________________
1796 Int_t AliTPCtrackerMI::FollowProlongationFast(AliTPCseed& t, Int_t rf, Int_t step) {
1797 //-----------------------------------------------------------------
1798 // This function tries to find a track prolongation.
1799 //-----------------------------------------------------------------
1800 Double_t xt=t.GetX();
1802 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
1803 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1804 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1805 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
1807 for (Int_t nr=GetRowNumber(xt)-1; nr>=rf; nr-=step) {
1809 if (FollowToNextFast(t,nr)==0)
1810 if (!t.IsActive()) return 0;
1820 Int_t AliTPCtrackerMI::FollowBackProlongation(AliTPCseed& t, Int_t rf) {
1821 //-----------------------------------------------------------------
1822 // This function tries to find a track prolongation.
1823 //-----------------------------------------------------------------
1825 Double_t xt=t.GetX();
1826 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
1827 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1828 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1829 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
1831 Int_t first = t.fFirstPoint;
1832 if (first<GetRowNumber(xt)+1) first = GetRowNumber(xt)+1;
1834 if (first<0) first=0;
1835 for (Int_t nr=first; nr<=rf; nr++) {
1836 if ( (TMath::Abs(t.GetSnp())>0.95)) break;
1837 if (t.GetKinkIndexes()[0]<0){
1838 for (Int_t i=0;i<3;i++){
1839 Int_t index = t.GetKinkIndexes()[i];
1840 if (index==0) break;
1841 if (index>0) continue;
1842 index = TMath::Abs(index);
1843 AliESDkink * kink = fEvent->GetKink(index-1);
1845 printf("PROBLEM\n");
1848 Int_t kinkrow = kink->GetTPCRow0()-2-Int_t(0.5/(0.05+kink->GetAngle(2)));
1850 AliExternalTrackParam paramm(t);
1851 kink->SetMother(paramm);
1852 kink->SetStatus(2,1);
1859 if (nr<fInnerSec->GetNRows())
1860 fSectors = fInnerSec;
1862 fSectors = fOuterSec;
1872 Float_t AliTPCtrackerMI::OverlapFactor(AliTPCseed * s1, AliTPCseed * s2, Int_t &sum1, Int_t & sum2)
1880 Float_t dz2 =(s1->GetZ() - s2->GetZ());
1883 Float_t dy2 =TMath::Abs((s1->GetY() - s2->GetY()));
1885 Float_t distance = TMath::Sqrt(dz2+dy2);
1886 if (distance>4.) return 0; // if there are far away - not overlap - to reduce combinatorics
1889 Int_t firstpoint = TMath::Min(s1->fFirstPoint,s2->fFirstPoint);
1890 Int_t lastpoint = TMath::Max(s1->fLastPoint,s2->fLastPoint);
1895 if (firstpoint>lastpoint) {
1896 firstpoint =lastpoint;
1901 for (Int_t i=firstpoint-1;i<lastpoint+1;i++){
1902 if (s1->GetClusterIndex2(i)>0) sum1++;
1903 if (s2->GetClusterIndex2(i)>0) sum2++;
1904 if (s1->GetClusterIndex2(i)==s2->GetClusterIndex2(i) && s1->GetClusterIndex2(i)>0) {
1908 if (sum<5) return 0;
1910 Float_t summin = TMath::Min(sum1+1,sum2+1);
1911 Float_t ratio = (sum+1)/Float_t(summin);
1915 void AliTPCtrackerMI::SignShared(AliTPCseed * s1, AliTPCseed * s2)
1919 if (TMath::Abs(s1->GetC()-s2->GetC())>0.004) return;
1920 if (TMath::Abs(s1->GetTgl()-s2->GetTgl())>0.6) return;
1922 Float_t dz2 =(s1->GetZ() - s2->GetZ());
1924 Float_t dy2 =(s1->GetY() - s2->GetY());
1926 Float_t distance = dz2+dy2;
1927 if (distance>325.) return ; // if there are far away - not overlap - to reduce combinatorics
1932 Int_t firstpoint = TMath::Max(s1->fFirstPoint,s2->fFirstPoint);
1933 Int_t lastpoint = TMath::Min(s1->fLastPoint,s2->fLastPoint);
1935 if (firstpoint>=lastpoint-5) return;;
1937 for (Int_t i=firstpoint;i<lastpoint;i++){
1938 // if ( (s1->GetClusterIndex2(i)&0xFFFF8FFF)==(s2->GetClusterIndex2(i)&0xFFFF8FFF) && s1->GetClusterIndex2(i)>0) {
1939 if ( (s1->GetClusterIndex2(i))==(s2->GetClusterIndex2(i)) && s1->GetClusterIndex2(i)>0) {
1946 for (Int_t i=firstpoint;i<lastpoint;i++){
1947 // if ( (s1->GetClusterIndex2(i)&0xFFFF8FFF)==(s2->GetClusterIndex2(i)&0xFFFF8FFF) && s1->GetClusterIndex2(i)>0) {
1948 if ( (s1->GetClusterIndex2(i))==(s2->GetClusterIndex2(i)) && s1->GetClusterIndex2(i)>0) {
1949 AliTPCTrackerPoint *p1 = s1->GetTrackPoint(i);
1950 AliTPCTrackerPoint *p2 = s2->GetTrackPoint(i);;
1951 if (s1->IsActive()&&s2->IsActive()){
1952 p1->fIsShared = kTRUE;
1953 p2->fIsShared = kTRUE;
1960 for (Int_t i=0;i<4;i++){
1961 if (s1->fOverlapLabels[3*i]==0){
1962 s1->fOverlapLabels[3*i] = s2->GetLabel();
1963 s1->fOverlapLabels[3*i+1] = sumshared;
1964 s1->fOverlapLabels[3*i+2] = s2->GetUniqueID();
1968 for (Int_t i=0;i<4;i++){
1969 if (s2->fOverlapLabels[3*i]==0){
1970 s2->fOverlapLabels[3*i] = s1->GetLabel();
1971 s2->fOverlapLabels[3*i+1] = sumshared;
1972 s2->fOverlapLabels[3*i+2] = s1->GetUniqueID();
1980 void AliTPCtrackerMI::SignShared(TObjArray * arr)
1983 //sort trackss according sectors
1985 for (Int_t i=0; i<arr->GetEntriesFast(); i++) {
1986 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
1988 //if (pt) RotateToLocal(pt);
1992 arr->Sort(); // sorting according z
1993 arr->Expand(arr->GetEntries());
1996 Int_t nseed=arr->GetEntriesFast();
1997 for (Int_t i=0; i<nseed; i++) {
1998 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2000 for (Int_t j=0;j<=12;j++){
2001 pt->fOverlapLabels[j] =0;
2004 for (Int_t i=0; i<nseed; i++) {
2005 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2007 if (pt->fRemoval>10) continue;
2008 for (Int_t j=i+1; j<nseed; j++){
2009 AliTPCseed *pt2=(AliTPCseed*)arr->UncheckedAt(j);
2011 if (pt2->fRemoval<=10) {
2012 if ( TMath::Abs(pt->fRelativeSector-pt2->fRelativeSector)>0) break;
2019 void AliTPCtrackerMI::RemoveDouble(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t removalindex)
2022 //sort trackss according sectors
2025 Info("RemoveDouble","Number of tracks before double removal- %d\n",arr->GetEntries());
2028 for (Int_t i=0; i<arr->GetEntriesFast(); i++) {
2029 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2034 arr->Sort(); // sorting according z
2035 arr->Expand(arr->GetEntries());
2037 //reset overlap labels
2039 Int_t nseed=arr->GetEntriesFast();
2040 for (Int_t i=0; i<nseed; i++) {
2041 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2044 for (Int_t j=0;j<=12;j++){
2045 pt->fOverlapLabels[j] =0;
2049 //sign shared tracks
2050 for (Int_t i=0; i<nseed; i++) {
2051 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2053 if (pt->fRemoval>10) continue;
2054 Float_t deltac = pt->GetC()*0.1;
2055 for (Int_t j=i+1; j<nseed; j++){
2056 AliTPCseed *pt2=(AliTPCseed*)arr->UncheckedAt(j);
2058 if (pt2->fRemoval<=10) {
2059 if ( TMath::Abs(pt->fRelativeSector-pt2->fRelativeSector)>0) break;
2060 if (TMath::Abs(pt->GetC() -pt2->GetC())>deltac) continue;
2061 if (TMath::Abs(pt->GetTgl()-pt2->GetTgl())>0.05) continue;
2068 // remove highly shared tracks
2069 for (Int_t i=0; i<nseed; i++) {
2070 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2072 if (pt->fRemoval>10) continue;
2075 for (Int_t j=0;j<4;j++){
2076 sumshared = pt->fOverlapLabels[j*3+1];
2078 Float_t factor = factor1;
2079 if (pt->fRemoval>0) factor = factor2;
2080 if (sumshared/pt->GetNumberOfClusters()>factor){
2081 for (Int_t j=0;j<4;j++){
2082 if (pt->fOverlapLabels[3*j]==0) continue;
2083 if (pt->fOverlapLabels[3*j+1]<5) continue;
2084 if (pt->fRemoval==removalindex) continue;
2085 AliTPCseed * pt2 = (AliTPCseed*)arr->UncheckedAt(pt->fOverlapLabels[3*j+2]);
2087 if (pt2->GetSigma2C()<pt->GetSigma2C()){
2088 // pt->fRemoval = removalindex;
2089 delete arr->RemoveAt(i);
2097 Info("RemoveDouble","Number of tracks after double removal- %d\n",arr->GetEntries());
2106 void AliTPCtrackerMI::SortTracks(TObjArray * arr, Int_t mode) const
2109 //sort tracks in array according mode criteria
2110 Int_t nseed = arr->GetEntriesFast();
2111 for (Int_t i=0; i<nseed; i++) {
2112 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2122 void AliTPCtrackerMI::RemoveUsed(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t removalindex)
2125 //Loop over all tracks and remove "overlaps"
2128 Int_t nseed = arr->GetEntriesFast();
2131 for (Int_t i=0; i<nseed; i++) {
2132 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2134 delete arr->RemoveAt(i);
2138 pt->fBSigned = kFALSE;
2142 nseed = arr->GetEntriesFast();
2149 for (Int_t i=0; i<nseed; i++) {
2150 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2154 Int_t found,foundable,shared;
2156 pt->GetClusterStatistic(0,160,found, foundable,shared,kFALSE);
2158 pt->GetClusterStatistic(0,160,found, foundable,shared,kTRUE);
2160 Double_t factor = factor2;
2161 if (pt->fBConstrain) factor = factor1;
2163 if ((Float_t(shared)/Float_t(found))>factor){
2164 pt->Desactivate(removalindex);
2169 for (Int_t i=0; i<160; i++) {
2170 Int_t index=pt->GetClusterIndex2(i);
2171 if (index<0 || index&0x8000 ) continue;
2172 AliTPCclusterMI *c= pt->fClusterPointer[i];
2174 // if (!c->IsUsed(10)) c->Use(10);
2175 //if (pt->IsActive())
2184 Info("RemoveUsed","\n*****\nNumber of good tracks after shared removal\t%d\n",fNtracks);
2189 void AliTPCtrackerMI::RemoveUsed2(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t minimal)
2192 //Loop over all tracks and remove "overlaps"
2197 Int_t nseed = arr->GetEntriesFast();
2198 Float_t * quality = new Float_t[nseed];
2199 Int_t * indexes = new Int_t[nseed];
2203 for (Int_t i=0; i<nseed; i++) {
2204 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2209 pt->UpdatePoints(); //select first last max dens points
2210 Float_t * points = pt->GetPoints();
2211 if (points[3]<0.8) quality[i] =-1;
2213 quality[i] = (points[2]-points[0])+pt->GetNumberOfClusters();
2215 TMath::Sort(nseed,quality,indexes);
2218 for (Int_t itrack=0; itrack<nseed; itrack++) {
2219 Int_t trackindex = indexes[itrack];
2220 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(trackindex);
2221 if (quality[trackindex]<0){
2223 delete arr->RemoveAt(trackindex);
2226 arr->RemoveAt(trackindex);
2231 Int_t first = Int_t(pt->GetPoints()[0]);
2232 Int_t last = Int_t(pt->GetPoints()[2]);
2233 Double_t factor = (pt->fBConstrain) ? factor1: factor2;
2235 Int_t found,foundable,shared;
2236 pt->GetClusterStatistic(first,last, found, foundable,shared,kFALSE);
2237 Float_t sharedfactor = Float_t(shared+1)/Float_t(found+1);
2238 Bool_t itsgold =kFALSE;
2241 if (pt->fEsd->GetITSclusters(dummy)>4) itsgold= kTRUE;
2245 if (Float_t(shared+1)/Float_t(found+1)>factor){
2246 if (pt->GetKinkIndexes()[0]!=0) continue; //don't remove tracks - part of the kinks
2247 delete arr->RemoveAt(trackindex);
2250 if (pt->GetNumberOfClusters()<50&&(found-0.5*shared)<minimal){ //remove short tracks
2251 if (pt->GetKinkIndexes()[0]!=0) continue; //don't remove tracks - part of the kinks
2252 delete arr->RemoveAt(trackindex);
2258 if (sharedfactor>0.4) continue;
2259 if (pt->GetKinkIndexes()[0]>0) continue;
2260 for (Int_t i=first; i<last; i++) {
2261 Int_t index=pt->GetClusterIndex2(i);
2262 // if (index<0 || index&0x8000 ) continue;
2263 if (index<0 || index&0x8000 ) continue;
2264 AliTPCclusterMI *c= pt->fClusterPointer[i];
2271 Info("RemoveUsed","\n*****\nNumber of good tracks after shared removal\t%d\n",fNtracks);
2277 void AliTPCtrackerMI::UnsignClusters()
2280 // loop over all clusters and unsign them
2283 for (Int_t sec=0;sec<fkNIS;sec++){
2284 for (Int_t row=0;row<fInnerSec->GetNRows();row++){
2285 AliTPCclusterMI *cl = fInnerSec[sec][row].fClusters1;
2286 for (Int_t icl =0;icl< fInnerSec[sec][row].fN1;icl++)
2287 // if (cl[icl].IsUsed(10))
2289 cl = fInnerSec[sec][row].fClusters2;
2290 for (Int_t icl =0;icl< fInnerSec[sec][row].fN2;icl++)
2291 //if (cl[icl].IsUsed(10))
2296 for (Int_t sec=0;sec<fkNOS;sec++){
2297 for (Int_t row=0;row<fOuterSec->GetNRows();row++){
2298 AliTPCclusterMI *cl = fOuterSec[sec][row].fClusters1;
2299 for (Int_t icl =0;icl< fOuterSec[sec][row].fN1;icl++)
2300 //if (cl[icl].IsUsed(10))
2302 cl = fOuterSec[sec][row].fClusters2;
2303 for (Int_t icl =0;icl< fOuterSec[sec][row].fN2;icl++)
2304 //if (cl[icl].IsUsed(10))
2313 void AliTPCtrackerMI::SignClusters(TObjArray * arr, Float_t fnumber, Float_t fdensity)
2316 //sign clusters to be "used"
2318 // snumber and sdensity sign number of sigmas - bellow mean value to be accepted
2319 // loop over "primaries"
2333 Int_t nseed = arr->GetEntriesFast();
2334 for (Int_t i=0; i<nseed; i++) {
2335 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2339 if (!(pt->IsActive())) continue;
2340 Float_t dens = pt->GetNumberOfClusters()/Float_t(pt->fNFoundable);
2341 if ( (dens>0.7) && (pt->GetNumberOfClusters()>70)){
2343 sumdens2+= dens*dens;
2344 sumn += pt->GetNumberOfClusters();
2345 sumn2 += pt->GetNumberOfClusters()*pt->GetNumberOfClusters();
2346 Float_t chi2 = pt->GetChi2()/pt->GetNumberOfClusters();
2349 sumchi2 +=chi2*chi2;
2354 Float_t mdensity = 0.9;
2355 Float_t meann = 130;
2356 Float_t meanchi = 1;
2357 Float_t sdensity = 0.1;
2358 Float_t smeann = 10;
2359 Float_t smeanchi =0.4;
2363 mdensity = sumdens/sum;
2365 meanchi = sumchi/sum;
2367 sdensity = sumdens2/sum-mdensity*mdensity;
2368 sdensity = TMath::Sqrt(sdensity);
2370 smeann = sumn2/sum-meann*meann;
2371 smeann = TMath::Sqrt(smeann);
2373 smeanchi = sumchi2/sum - meanchi*meanchi;
2374 smeanchi = TMath::Sqrt(smeanchi);
2378 //REMOVE SHORT DELTAS or tracks going out of sensitive volume of TPC
2380 for (Int_t i=0; i<nseed; i++) {
2381 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2385 if (pt->fBSigned) continue;
2386 if (pt->fBConstrain) continue;
2387 //if (!(pt->IsActive())) continue;
2389 Int_t found,foundable,shared;
2390 pt->GetClusterStatistic(0,160,found, foundable,shared);
2391 if (shared/float(found)>0.3) {
2392 if (shared/float(found)>0.9 ){
2393 //delete arr->RemoveAt(i);
2398 Bool_t isok =kFALSE;
2399 if ( (pt->fNShared/pt->GetNumberOfClusters()<0.5) &&pt->GetNumberOfClusters()>60)
2401 if ((TMath::Abs(1/pt->GetC())<100.) && (pt->fNShared/pt->GetNumberOfClusters()<0.7))
2403 if (TMath::Abs(pt->GetZ()/pt->GetX())>1.1)
2405 if ( (TMath::Abs(pt->GetSnp()>0.7) && pt->GetD(0,0)>60.))
2409 for (Int_t i=0; i<160; i++) {
2410 Int_t index=pt->GetClusterIndex2(i);
2411 if (index<0) continue;
2412 AliTPCclusterMI *c= pt->fClusterPointer[i];
2414 //if (!(c->IsUsed(10))) c->Use();
2421 Double_t maxchi = meanchi+2.*smeanchi;
2423 for (Int_t i=0; i<nseed; i++) {
2424 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2428 //if (!(pt->IsActive())) continue;
2429 if (pt->fBSigned) continue;
2430 Double_t chi = pt->GetChi2()/pt->GetNumberOfClusters();
2431 if (chi>maxchi) continue;
2434 Float_t dens = pt->GetNumberOfClusters()/Float_t(pt->fNFoundable);
2436 //sign only tracks with enoug big density at the beginning
2438 if ((pt->GetDensityFirst(40)<0.75) && pt->GetNumberOfClusters()<meann) continue;
2441 Double_t mindens = TMath::Max(double(mdensity-sdensity*fdensity*bfactor),0.65);
2442 Double_t minn = TMath::Max(Int_t(meann-fnumber*smeann*bfactor),50);
2444 // if (pt->fBConstrain) mindens = TMath::Max(mdensity-sdensity*fdensity*bfactor,0.65);
2445 if ( (pt->fRemoval==10) && (pt->GetSnp()>0.8)&&(dens>mindens))
2448 if ((dens>mindens && pt->GetNumberOfClusters()>minn) && chi<maxchi ){
2449 //Int_t noc=pt->GetNumberOfClusters();
2450 pt->fBSigned = kTRUE;
2451 for (Int_t i=0; i<160; i++) {
2453 Int_t index=pt->GetClusterIndex2(i);
2454 if (index<0) continue;
2455 AliTPCclusterMI *c= pt->fClusterPointer[i];
2457 // if (!(c->IsUsed(10))) c->Use();
2462 // gLastCheck = nseed;
2470 void AliTPCtrackerMI::StopNotActive(TObjArray * arr, Int_t row0, Float_t th0, Float_t th1, Float_t th2) const
2472 // stop not active tracks
2473 // take th1 as threshold for number of founded to number of foundable on last 10 active rows
2474 // take th2 as threshold for number of founded to number of foundable on last 20 active rows
2475 Int_t nseed = arr->GetEntriesFast();
2477 for (Int_t i=0; i<nseed; i++) {
2478 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2482 if (!(pt->IsActive())) continue;
2483 StopNotActive(pt,row0,th0, th1,th2);
2489 void AliTPCtrackerMI::StopNotActive(AliTPCseed * seed, Int_t row0, Float_t th0, Float_t th1,
2492 // stop not active tracks
2493 // take th1 as threshold for number of founded to number of foundable on last 10 active rows
2494 // take th2 as threshold for number of founded to number of foundable on last 20 active rows
2497 Int_t foundable = 0;
2498 Int_t maxindex = seed->fLastPoint; //last foundable row
2499 if (seed->fNFoundable*th0 > seed->GetNumberOfClusters()) {
2500 seed->Desactivate(10) ;
2504 for (Int_t i=row0; i<maxindex; i++){
2505 Int_t index = seed->GetClusterIndex2(i);
2506 if (index!=-1) foundable++;
2508 if (foundable<=30) sumgood1++;
2509 if (foundable<=50) {
2516 if (foundable>=30.){
2517 if (sumgood1<(th1*30.)) seed->Desactivate(10);
2520 if (sumgood2<(th2*50.)) seed->Desactivate(10);
2524 Int_t AliTPCtrackerMI::RefitInward(AliESD *event)
2527 // back propagation of ESD tracks
2533 //PrepareForProlongation(fSeeds,1);
2534 PropagateForward2(fSeeds);
2537 Int_t nseed = fSeeds->GetEntriesFast();
2538 for (Int_t i=0;i<nseed;i++){
2539 AliTPCseed * seed = (AliTPCseed*) fSeeds->UncheckedAt(i);
2540 if (!seed) continue;
2541 if (seed->GetKinkIndex(0)>0) UpdateKinkQualityD(seed); // update quality informations for kinks
2543 seed->PropagateTo(fParam->GetInnerRadiusLow());
2544 seed->UpdatePoints();
2545 AliESDtrack *esd=event->GetTrack(i);
2546 seed->CookdEdx(0.02,0.6);
2547 CookLabel(seed,0.1); //For comparison only
2549 if (0 && seed!=0&&esd!=0) {
2550 TTreeSRedirector &cstream = *fDebugStreamer;
2556 if (seed->GetNumberOfClusters()>15){
2557 esd->UpdateTrackParams(seed,AliESDtrack::kTPCrefit);
2558 esd->SetTPCPoints(seed->GetPoints());
2559 esd->SetTPCPointsF(seed->fNFoundable);
2560 Int_t ndedx = seed->fNCDEDX[0]+seed->fNCDEDX[1]+seed->fNCDEDX[2]+seed->fNCDEDX[3];
2561 Float_t sdedx = (seed->fSDEDX[0]+seed->fSDEDX[1]+seed->fSDEDX[2]+seed->fSDEDX[3])*0.25;
2562 Float_t dedx = seed->GetdEdx();
2563 esd->SetTPCsignal(dedx, sdedx, ndedx);
2567 //printf("problem\n");
2570 //FindKinks(fSeeds,event);
2571 Info("RefitInward","Number of refitted tracks %d",ntracks);
2578 Int_t AliTPCtrackerMI::PropagateBack(AliESD *event)
2581 // back propagation of ESD tracks
2587 PropagateBack(fSeeds);
2588 RemoveUsed2(fSeeds,0.4,0.4,20);
2590 Int_t nseed = fSeeds->GetEntriesFast();
2592 for (Int_t i=0;i<nseed;i++){
2593 AliTPCseed * seed = (AliTPCseed*) fSeeds->UncheckedAt(i);
2594 if (!seed) continue;
2595 if (seed->GetKinkIndex(0)<0) UpdateKinkQualityM(seed); // update quality informations for kinks
2596 seed->UpdatePoints();
2597 AliESDtrack *esd=event->GetTrack(i);
2598 seed->CookdEdx(0.02,0.6);
2599 CookLabel(seed,0.1); //For comparison only
2600 if (seed->GetNumberOfClusters()>15){
2601 esd->UpdateTrackParams(seed,AliESDtrack::kTPCout);
2602 esd->SetTPCPoints(seed->GetPoints());
2606 //FindKinks(fSeeds,event);
2607 Info("PropagateBack","Number of back propagated tracks %d",ntracks);
2614 void AliTPCtrackerMI::DeleteSeeds()
2618 Int_t nseed = fSeeds->GetEntriesFast();
2619 for (Int_t i=0;i<nseed;i++){
2620 AliTPCseed * seed = (AliTPCseed*)fSeeds->At(i);
2621 if (seed) delete fSeeds->RemoveAt(i);
2627 void AliTPCtrackerMI::ReadSeeds(AliESD *event, Int_t direction)
2630 //read seeds from the event
2632 Int_t nentr=event->GetNumberOfTracks();
2634 Info("ReadSeeds", "Number of ESD tracks: %d\n", nentr);
2639 fSeeds = new TObjArray(nentr);
2643 for (Int_t i=0; i<nentr; i++) {
2644 AliESDtrack *esd=event->GetTrack(i);
2645 ULong_t status=esd->GetStatus();
2646 if (!(status&AliESDtrack::kTPCin)) continue;
2647 AliTPCtrack t(*esd);
2648 t.SetNumberOfClusters(0);
2649 // AliTPCseed *seed = new AliTPCseed(t,t.GetAlpha());
2650 AliTPCseed *seed = new AliTPCseed(t/*,t.GetAlpha()*/);
2651 for (Int_t ikink=0;ikink<3;ikink++) {
2652 Int_t index = esd->GetKinkIndex(ikink);
2653 seed->GetKinkIndexes()[ikink] = index;
2654 if (index==0) continue;
2655 index = TMath::Abs(index);
2656 AliESDkink * kink = fEvent->GetKink(index-1);
2657 if (kink&&esd->GetKinkIndex(ikink)<0){
2658 if ((status & AliESDtrack::kTRDrefit) != 0) kink->SetStatus(1,2);
2659 if ((status & AliESDtrack::kITSout) != 0) kink->SetStatus(1,0);
2661 if (kink&&esd->GetKinkIndex(ikink)>0){
2662 if ((status & AliESDtrack::kTRDrefit) != 0) kink->SetStatus(1,6);
2663 if ((status & AliESDtrack::kITSout) != 0) kink->SetStatus(1,4);
2667 if (((status&AliESDtrack::kITSout)==0)&&(direction==1)) seed->ResetCovariance();
2668 if ( direction ==2 &&(status & AliESDtrack::kTRDrefit) == 0 ) seed->ResetCovariance();
2669 if ( direction ==2 && ((status & AliESDtrack::kTPCout) == 0) ) {
2674 if ( direction ==2 &&(status & AliESDtrack::kTRDrefit) > 0 ) {
2675 Double_t par0[5],par1[5],alpha,x;
2676 esd->GetInnerExternalParameters(alpha,x,par0);
2677 esd->GetExternalParameters(x,par1);
2678 Double_t delta1 = TMath::Abs(par0[4]-par1[4])/(0.000000001+TMath::Abs(par0[4]+par1[4]));
2679 Double_t delta2 = TMath::Abs(par0[3]-par1[3]);
2681 if (esd->GetTRDncls()>0) trdchi2 = esd->GetTRDchi2()/esd->GetTRDncls();
2682 //reset covariance if suspicious
2683 if ( (delta1>0.1) || (delta2>0.006) ||trdchi2>7.)
2684 seed->ResetCovariance();
2689 // rotate to the local coordinate system
2691 fSectors=fInnerSec; fN=fkNIS;
2692 Double_t alpha=seed->GetAlpha() - fSectors->GetAlphaShift();
2693 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
2694 if (alpha < 0. ) alpha += 2.*TMath::Pi();
2695 Int_t ns=Int_t(alpha/fSectors->GetAlpha())%fN;
2696 alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift();
2697 if (alpha<-TMath::Pi()) alpha += 2*TMath::Pi();
2698 if (alpha>=TMath::Pi()) alpha -= 2*TMath::Pi();
2699 alpha-=seed->GetAlpha();
2700 if (!seed->Rotate(alpha)) {
2706 if (esd->GetKinkIndex(0)<=0){
2707 for (Int_t irow=0;irow<160;irow++){
2708 Int_t index = seed->GetClusterIndex2(irow);
2711 AliTPCclusterMI * cl = GetClusterMI(index);
2712 seed->fClusterPointer[irow] = cl;
2714 if ((index & 0x8000)==0){
2715 cl->Use(10); // accepted cluster
2717 cl->Use(6); // close cluster not accepted
2720 Info("ReadSeeds","Not found cluster");
2725 fSeeds->AddAt(seed,i);
2731 //_____________________________________________________________________________
2732 void AliTPCtrackerMI::MakeSeeds3(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4],
2733 Float_t deltay, Int_t ddsec) {
2734 //-----------------------------------------------------------------
2735 // This function creates track seeds.
2736 // SEEDING WITH VERTEX CONSTRAIN
2737 //-----------------------------------------------------------------
2738 // cuts[0] - fP4 cut
2739 // cuts[1] - tan(phi) cut
2740 // cuts[2] - zvertex cut
2741 // cuts[3] - fP3 cut
2749 Double_t x[5], c[15];
2750 // Int_t di = i1-i2;
2752 AliTPCseed * seed = new AliTPCseed;
2753 Double_t alpha=fSectors->GetAlpha(), shift=fSectors->GetAlphaShift();
2754 Double_t cs=cos(alpha), sn=sin(alpha);
2756 // Double_t x1 =fOuterSec->GetX(i1);
2757 //Double_t xx2=fOuterSec->GetX(i2);
2759 Double_t x1 =GetXrow(i1);
2760 Double_t xx2=GetXrow(i2);
2762 Double_t x3=GetX(), y3=GetY(), z3=GetZ();
2764 Int_t imiddle = (i2+i1)/2; //middle pad row index
2765 Double_t xm = GetXrow(imiddle); // radius of middle pad-row
2766 const AliTPCRow& krm=GetRow(sec,imiddle); //middle pad -row
2770 const AliTPCRow& kr1=GetRow(ns,i1);
2771 Double_t ymax = GetMaxY(i1)-kr1.fDeadZone-1.5;
2772 Double_t ymaxm = GetMaxY(imiddle)-kr1.fDeadZone-1.5;
2775 // change cut on curvature if it can't reach this layer
2776 // maximal curvature set to reach it
2777 Double_t dvertexmax = TMath::Sqrt((x1-x3)*(x1-x3)+(ymax+5-y3)*(ymax+5-y3));
2778 if (dvertexmax*0.5*cuts[0]>0.85){
2779 cuts[0] = 0.85/(dvertexmax*0.5+1.);
2781 Double_t r2min = 1/(cuts[0]*cuts[0]); //minimal square of radius given by cut
2784 if (deltay>0) ddsec = 0;
2785 // loop over clusters
2786 for (Int_t is=0; is < kr1; is++) {
2788 if (kr1[is]->IsUsed(10)) continue;
2789 Double_t y1=kr1[is]->GetY(), z1=kr1[is]->GetZ();
2790 //if (TMath::Abs(y1)>ymax) continue;
2792 if (deltay>0 && TMath::Abs(ymax-TMath::Abs(y1))> deltay ) continue; // seed only at the edge
2794 // find possible directions
2795 Float_t anglez = (z1-z3)/(x1-x3);
2796 Float_t extraz = z1 - anglez*(x1-xx2); // extrapolated z
2799 //find rotation angles relative to line given by vertex and point 1
2800 Double_t dvertex2 = (x1-x3)*(x1-x3)+(y1-y3)*(y1-y3);
2801 Double_t dvertex = TMath::Sqrt(dvertex2);
2802 Double_t angle13 = TMath::ATan((y1-y3)/(x1-x3));
2803 Double_t cs13 = cos(-angle13), sn13 = sin(-angle13);
2806 // loop over 2 sectors
2812 Double_t dddz1=0; // direction of delta inclination in z axis
2819 for (Int_t dsec = dsec1; dsec<=dsec2;dsec++){
2820 Int_t sec2 = sec + dsec;
2822 // AliTPCRow& kr2 = fOuterSec[(sec2+fkNOS)%fkNOS][i2];
2823 //AliTPCRow& kr2m = fOuterSec[(sec2+fkNOS)%fkNOS][imiddle];
2824 AliTPCRow& kr2 = GetRow((sec2+fkNOS)%fkNOS,i2);
2825 AliTPCRow& kr2m = GetRow((sec2+fkNOS)%fkNOS,imiddle);
2826 Int_t index1 = TMath::Max(kr2.Find(extraz-0.6-dddz1*TMath::Abs(z1)*0.05)-1,0);
2827 Int_t index2 = TMath::Min(kr2.Find(extraz+0.6+dddz2*TMath::Abs(z1)*0.05)+1,kr2);
2829 // rotation angles to p1-p3
2830 Double_t cs13r = cos(-angle13+dsec*alpha)/dvertex, sn13r = sin(-angle13+dsec*alpha)/dvertex;
2831 Double_t x2, y2, z2;
2833 // Double_t dymax = maxangle*TMath::Abs(x1-xx2);
2836 Double_t dxx0 = (xx2-x3)*cs13r;
2837 Double_t dyy0 = (xx2-x3)*sn13r;
2838 for (Int_t js=index1; js < index2; js++) {
2839 const AliTPCclusterMI *kcl = kr2[js];
2840 if (kcl->IsUsed(10)) continue;
2842 //calcutate parameters
2844 Double_t yy0 = dyy0 +(kcl->GetY()-y3)*cs13r;
2846 if (TMath::Abs(yy0)<0.000001) continue;
2847 Double_t xx0 = dxx0 -(kcl->GetY()-y3)*sn13r;
2848 Double_t y0 = 0.5*(xx0*xx0+yy0*yy0-xx0)/yy0;
2849 Double_t r02 = (0.25+y0*y0)*dvertex2;
2850 //curvature (radius) cut
2851 if (r02<r2min) continue;
2855 Double_t c0 = 1/TMath::Sqrt(r02);
2859 //Double_t dfi0 = 2.*TMath::ASin(dvertex*c0*0.5);
2860 //Double_t dfi1 = 2.*TMath::ASin(TMath::Sqrt(yy0*yy0+(1-xx0)*(1-xx0))*dvertex*c0*0.5);
2861 Double_t dfi0 = 2.*AliTPCFastMath::FastAsin(dvertex*c0*0.5);
2862 Double_t dfi1 = 2.*AliTPCFastMath::FastAsin(TMath::Sqrt(yy0*yy0+(1-xx0)*(1-xx0))*dvertex*c0*0.5);
2865 Double_t z0 = kcl->GetZ();
2866 Double_t zzzz2 = z1-(z1-z3)*dfi1/dfi0;
2867 if (TMath::Abs(zzzz2-z0)>0.5) continue;
2870 Double_t dip = (z1-z0)*c0/dfi1;
2871 Double_t x0 = (0.5*cs13+y0*sn13)*dvertex*c0;
2882 x2= xx2*cs-y2*sn*dsec;
2883 y2=+xx2*sn*dsec+y2*cs;
2893 // do we have cluster at the middle ?
2895 GetProlongation(x1,xm,x,ym,zm);
2897 AliTPCclusterMI * cm=0;
2898 if (TMath::Abs(ym)-ymaxm<0){
2899 cm = krm.FindNearest2(ym,zm,1.0,0.6,dummy);
2900 if ((!cm) || (cm->IsUsed(10))) {
2905 // rotate y1 to system 0
2906 // get state vector in rotated system
2907 Double_t yr1 = (-0.5*sn13+y0*cs13)*dvertex*c0;
2908 Double_t xr2 = x0*cs+yr1*sn*dsec;
2909 Double_t xr[5]={kcl->GetY(),kcl->GetZ(), xr2, dip, c0};
2911 GetProlongation(xx2,xm,xr,ym,zm);
2912 if (TMath::Abs(ym)-ymaxm<0){
2913 cm = kr2m.FindNearest2(ym,zm,1.0,0.6,dummy);
2914 if ((!cm) || (cm->IsUsed(10))) {
2924 dym = ym - cm->GetY();
2925 dzm = zm - cm->GetZ();
2932 Double_t sy1=kr1[is]->GetSigmaY2()*2., sz1=kr1[is]->GetSigmaZ2()*2.;
2933 Double_t sy2=kcl->GetSigmaY2()*2., sz2=kcl->GetSigmaZ2()*2.;
2934 //Double_t sy3=400*3./12., sy=0.1, sz=0.1;
2935 Double_t sy3=25000*x[4]*x[4]+0.1, sy=0.1, sz=0.1;
2936 //Double_t sy3=25000*x[4]*x[4]*60+0.5, sy=0.1, sz=0.1;
2938 Double_t f40=(F1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
2939 Double_t f42=(F1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
2940 Double_t f43=(F1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
2941 Double_t f20=(F2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
2942 Double_t f22=(F2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
2943 Double_t f23=(F2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
2945 Double_t f30=(F3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
2946 Double_t f31=(F3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
2947 Double_t f32=(F3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
2948 Double_t f34=(F3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
2952 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
2953 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
2954 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
2955 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
2956 c[13]=f30*sy1*f40+f32*sy2*f42;
2957 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
2959 // if (!BuildSeed(kr1[is],kcl,0,x1,x2,x3,x,c)) continue;
2961 UInt_t index=kr1.GetIndex(is);
2962 AliTPCseed *track=new(seed) AliTPCseed(index, x, c, x1, ns*alpha+shift);
2964 track->fIsSeeding = kTRUE;
2971 FollowProlongation(*track, (i1+i2)/2,1);
2972 Int_t foundable,found,shared;
2973 track->GetClusterStatistic((i1+i2)/2,i1, found, foundable, shared, kTRUE);
2974 if ((found<0.55*foundable) || shared>0.5*found || (track->GetSigmaY2()+track->GetSigmaZ2())>0.5){
2976 seed->~AliTPCseed();
2982 FollowProlongation(*track, i2,1);
2986 track->fBConstrain =1;
2987 // track->fLastPoint = i1+fInnerSec->GetNRows(); // first cluster in track position
2988 track->fLastPoint = i1; // first cluster in track position
2989 track->fFirstPoint = track->fLastPoint;
2991 if (track->GetNumberOfClusters()<(i1-i2)*0.5 ||
2992 track->GetNumberOfClusters() < track->fNFoundable*0.6 ||
2993 track->fNShared>0.4*track->GetNumberOfClusters() ) {
2995 seed->~AliTPCseed();
2999 // Z VERTEX CONDITION
3001 zv = track->GetZ()+track->GetTgl()/track->GetC()*
3002 ( asin(-track->GetEta()) - asin(track->GetX()*track->GetC()-track->GetEta()));
3003 if (TMath::Abs(zv-z3)>cuts[2]) {
3004 FollowProlongation(*track, TMath::Max(i2-20,0));
3005 zv = track->GetZ()+track->GetTgl()/track->GetC()*
3006 ( asin(-track->GetEta()) - asin(track->GetX()*track->GetC()-track->GetEta()));
3007 if (TMath::Abs(zv-z3)>cuts[2]){
3008 FollowProlongation(*track, TMath::Max(i2-40,0));
3009 zv = track->GetZ()+track->GetTgl()/track->GetC()*
3010 ( asin(-track->GetEta()) - asin(track->GetX()*track->GetC()-track->GetEta()));
3011 if (TMath::Abs(zv-z3)>cuts[2] &&(track->GetNumberOfClusters() > track->fNFoundable*0.7)){
3012 // make seed without constrain
3013 AliTPCseed * track2 = MakeSeed(track,0.2,0.5,1.);
3014 FollowProlongation(*track2, i2,1);
3015 track2->fBConstrain = kFALSE;
3016 track2->fSeedType = 1;
3017 arr->AddLast(track2);
3019 seed->~AliTPCseed();
3024 seed->~AliTPCseed();
3031 track->fSeedType =0;
3032 arr->AddLast(track);
3033 seed = new AliTPCseed;
3035 // don't consider other combinations
3036 if (track->GetNumberOfClusters() > track->fNFoundable*0.8)
3042 Info("MakeSeeds3","\nSeeding statistic:\t%d\t%d\t%d\t%d\t%d\t%d",nin0,nin1,nin2,nin,nout1,nout2);
3048 void AliTPCtrackerMI::MakeSeeds5(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4],
3053 //-----------------------------------------------------------------
3054 // This function creates track seeds.
3055 //-----------------------------------------------------------------
3056 // cuts[0] - fP4 cut
3057 // cuts[1] - tan(phi) cut
3058 // cuts[2] - zvertex cut
3059 // cuts[3] - fP3 cut
3069 Double_t x[5], c[15];
3071 // make temporary seed
3072 AliTPCseed * seed = new AliTPCseed;
3073 Double_t alpha=fOuterSec->GetAlpha(), shift=fOuterSec->GetAlphaShift();
3074 // Double_t cs=cos(alpha), sn=sin(alpha);
3079 Double_t x1 = GetXrow(i1-1);
3080 const AliTPCRow& kr1=GetRow(sec,i1-1);
3081 Double_t y1max = GetMaxY(i1-1)-kr1.fDeadZone-1.5;
3083 Double_t x1p = GetXrow(i1);
3084 const AliTPCRow& kr1p=GetRow(sec,i1);
3086 Double_t x1m = GetXrow(i1-2);
3087 const AliTPCRow& kr1m=GetRow(sec,i1-2);
3090 //last 3 padrow for seeding
3091 AliTPCRow& kr3 = GetRow((sec+fkNOS)%fkNOS,i1-7);
3092 Double_t x3 = GetXrow(i1-7);
3093 // Double_t y3max= GetMaxY(i1-7)-kr3.fDeadZone-1.5;
3095 AliTPCRow& kr3p = GetRow((sec+fkNOS)%fkNOS,i1-6);
3096 Double_t x3p = GetXrow(i1-6);
3098 AliTPCRow& kr3m = GetRow((sec+fkNOS)%fkNOS,i1-8);
3099 Double_t x3m = GetXrow(i1-8);
3104 Int_t im = i1-4; //middle pad row index
3105 Double_t xm = GetXrow(im); // radius of middle pad-row
3106 const AliTPCRow& krm=GetRow(sec,im); //middle pad -row
3107 // Double_t ymmax = GetMaxY(im)-kr1.fDeadZone-1.5;
3110 Double_t deltax = x1-x3;
3111 Double_t dymax = deltax*cuts[1];
3112 Double_t dzmax = deltax*cuts[3];
3114 // loop over clusters
3115 for (Int_t is=0; is < kr1; is++) {
3117 if (kr1[is]->IsUsed(10)) continue;
3118 Double_t y1=kr1[is]->GetY(), z1=kr1[is]->GetZ();
3120 if (deltay>0 && TMath::Abs(y1max-TMath::Abs(y1))> deltay ) continue; // seed only at the edge
3122 Int_t index1 = TMath::Max(kr3.Find(z1-dzmax)-1,0);
3123 Int_t index2 = TMath::Min(kr3.Find(z1+dzmax)+1,kr3);
3129 for (Int_t js=index1; js < index2; js++) {
3130 const AliTPCclusterMI *kcl = kr3[js];
3131 if (kcl->IsUsed(10)) continue;
3133 // apply angular cuts
3134 if (TMath::Abs(y1-y3)>dymax) continue;
3137 if (TMath::Abs(z1-z3)>dzmax) continue;
3139 Double_t angley = (y1-y3)/(x1-x3);
3140 Double_t anglez = (z1-z3)/(x1-x3);
3142 Double_t erry = TMath::Abs(angley)*(x1-x1m)*0.5+0.5;
3143 Double_t errz = TMath::Abs(anglez)*(x1-x1m)*0.5+0.5;
3145 Double_t yyym = angley*(xm-x1)+y1;
3146 Double_t zzzm = anglez*(xm-x1)+z1;
3148 const AliTPCclusterMI *kcm = krm.FindNearest2(yyym,zzzm,erry,errz,index);
3150 if (kcm->IsUsed(10)) continue;
3152 erry = TMath::Abs(angley)*(x1-x1m)*0.4+0.5;
3153 errz = TMath::Abs(anglez)*(x1-x1m)*0.4+0.5;
3160 // look around first
3161 const AliTPCclusterMI *kc1m = kr1m.FindNearest2(angley*(x1m-x1)+y1,
3167 if (kc1m->IsUsed(10)) used++;
3169 const AliTPCclusterMI *kc1p = kr1p.FindNearest2(angley*(x1p-x1)+y1,
3175 if (kc1p->IsUsed(10)) used++;
3177 if (used>1) continue;
3178 if (found<1) continue;
3182 const AliTPCclusterMI *kc3m = kr3m.FindNearest2(angley*(x3m-x3)+y3,
3188 if (kc3m->IsUsed(10)) used++;
3192 const AliTPCclusterMI *kc3p = kr3p.FindNearest2(angley*(x3p-x3)+y3,
3198 if (kc3p->IsUsed(10)) used++;
3202 if (used>1) continue;
3203 if (found<3) continue;
3213 x[4]=F1(x1,y1,x2,y2,x3,y3);
3214 //if (TMath::Abs(x[4]) >= cuts[0]) continue;
3217 x[2]=F2(x1,y1,x2,y2,x3,y3);
3220 x[3]=F3n(x1,y1,x2,y2,z1,z2,x[4]);
3221 //if (TMath::Abs(x[3]) > cuts[3]) continue;
3225 Double_t sy1=0.1, sz1=0.1;
3226 Double_t sy2=0.1, sz2=0.1;
3227 Double_t sy3=0.1, sy=0.1, sz=0.1;
3229 Double_t f40=(F1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
3230 Double_t f42=(F1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
3231 Double_t f43=(F1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
3232 Double_t f20=(F2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
3233 Double_t f22=(F2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
3234 Double_t f23=(F2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
3236 Double_t f30=(F3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
3237 Double_t f31=(F3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
3238 Double_t f32=(F3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
3239 Double_t f34=(F3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
3243 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
3244 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
3245 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
3246 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
3247 c[13]=f30*sy1*f40+f32*sy2*f42;
3248 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
3250 // if (!BuildSeed(kr1[is],kcl,0,x1,x2,x3,x,c)) continue;
3252 UInt_t index=kr1.GetIndex(is);
3253 AliTPCseed *track=new(seed) AliTPCseed(index, x, c, x1, sec*alpha+shift);
3255 track->fIsSeeding = kTRUE;
3258 FollowProlongation(*track, i1-7,1);
3259 if (track->GetNumberOfClusters() < track->fNFoundable*0.75 ||
3260 track->fNShared>0.6*track->GetNumberOfClusters() || ( track->GetSigmaY2()+ track->GetSigmaZ2())>0.6){
3262 seed->~AliTPCseed();
3268 FollowProlongation(*track, i2,1);
3269 track->fBConstrain =0;
3270 track->fLastPoint = i1+fInnerSec->GetNRows(); // first cluster in track position
3271 track->fFirstPoint = track->fLastPoint;
3273 if (track->GetNumberOfClusters()<(i1-i2)*0.5 ||
3274 track->GetNumberOfClusters()<track->fNFoundable*0.7 ||
3275 track->fNShared>2. || track->GetChi2()/track->GetNumberOfClusters()>6 || ( track->GetSigmaY2()+ track->GetSigmaZ2())>0.5 ) {
3277 seed->~AliTPCseed();
3282 FollowProlongation(*track, TMath::Max(i2-10,0),1);
3283 AliTPCseed * track2 = MakeSeed(track,0.2,0.5,0.9);
3284 FollowProlongation(*track2, i2,1);
3285 track2->fBConstrain = kFALSE;
3286 track2->fSeedType = 4;
3287 arr->AddLast(track2);
3289 seed->~AliTPCseed();
3293 //arr->AddLast(track);
3294 //seed = new AliTPCseed;
3300 Info("MakeSeeds5","\nSeeding statiistic:\t%d\t%d\t%d\t%d\t%d\t%d",nin0,nin1,nin2,nin,nout1,nout2,nout3);
3306 //_____________________________________________________________________________
3307 void AliTPCtrackerMI::MakeSeeds2(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t */*cuts[4]*/,
3308 Float_t deltay, Bool_t /*bconstrain*/) {
3309 //-----------------------------------------------------------------
3310 // This function creates track seeds - without vertex constraint
3311 //-----------------------------------------------------------------
3312 // cuts[0] - fP4 cut - not applied
3313 // cuts[1] - tan(phi) cut
3314 // cuts[2] - zvertex cut - not applied
3315 // cuts[3] - fP3 cut
3325 Double_t alpha=fOuterSec->GetAlpha(), shift=fOuterSec->GetAlphaShift();
3326 // Double_t cs=cos(alpha), sn=sin(alpha);
3327 Int_t row0 = (i1+i2)/2;
3328 Int_t drow = (i1-i2)/2;
3329 const AliTPCRow& kr0=fSectors[sec][row0];
3332 AliTPCpolyTrack polytrack;
3333 Int_t nclusters=fSectors[sec][row0];
3334 AliTPCseed * seed = new AliTPCseed;
3339 for (Int_t is=0; is < nclusters; is++) { //LOOP over clusters
3341 Int_t nfoundable =0;
3342 for (Int_t iter =1; iter<2; iter++){ //iterations
3343 const AliTPCRow& krm=fSectors[sec][row0-iter];
3344 const AliTPCRow& krp=fSectors[sec][row0+iter];
3345 const AliTPCclusterMI * cl= kr0[is];
3347 if (cl->IsUsed(10)) {
3353 Double_t x = kr0.GetX();
3354 // Initialization of the polytrack
3359 Double_t y0= cl->GetY();
3360 Double_t z0= cl->GetZ();
3364 Double_t ymax = fSectors->GetMaxY(row0)-kr0.fDeadZone-1.5;
3365 if (deltay>0 && TMath::Abs(ymax-TMath::Abs(y0))> deltay ) continue; // seed only at the edge
3367 erry = (0.5)*cl->GetSigmaY2()/TMath::Sqrt(cl->GetQ())*6;
3368 errz = (0.5)*cl->GetSigmaZ2()/TMath::Sqrt(cl->GetQ())*6;
3369 polytrack.AddPoint(x,y0,z0,erry, errz);
3372 if (cl->IsUsed(10)) sumused++;
3375 Float_t roady = (5*TMath::Sqrt(cl->GetSigmaY2()+0.2)+1.)*iter;
3376 Float_t roadz = (5*TMath::Sqrt(cl->GetSigmaZ2()+0.2)+1.)*iter;
3379 AliTPCclusterMI * cl1 = krm.FindNearest(y0,z0,roady,roadz);
3380 if (cl1 && TMath::Abs(ymax-TMath::Abs(y0))) {
3381 erry = (0.5)*cl1->GetSigmaY2()/TMath::Sqrt(cl1->GetQ())*3;
3382 errz = (0.5)*cl1->GetSigmaZ2()/TMath::Sqrt(cl1->GetQ())*3;
3383 if (cl1->IsUsed(10)) sumused++;
3384 polytrack.AddPoint(x,cl1->GetY(),cl1->GetZ(),erry,errz);
3388 AliTPCclusterMI * cl2 = krp.FindNearest(y0,z0,roady,roadz);
3390 erry = (0.5)*cl2->GetSigmaY2()/TMath::Sqrt(cl2->GetQ())*3;
3391 errz = (0.5)*cl2->GetSigmaZ2()/TMath::Sqrt(cl2->GetQ())*3;
3392 if (cl2->IsUsed(10)) sumused++;
3393 polytrack.AddPoint(x,cl2->GetY(),cl2->GetZ(),erry,errz);
3396 if (sumused>0) continue;
3398 polytrack.UpdateParameters();
3404 nfoundable = polytrack.GetN();
3405 nfound = nfoundable;
3407 for (Int_t ddrow = iter+1; ddrow<drow;ddrow++){
3408 Float_t maxdist = 0.8*(1.+3./(ddrow));
3409 for (Int_t delta = -1;delta<=1;delta+=2){
3410 Int_t row = row0+ddrow*delta;
3411 kr = &(fSectors[sec][row]);
3412 Double_t xn = kr->GetX();
3413 Double_t ymax = fSectors->GetMaxY(row)-kr->fDeadZone-1.5;
3414 polytrack.GetFitPoint(xn,yn,zn);
3415 if (TMath::Abs(yn)>ymax) continue;
3417 AliTPCclusterMI * cln = kr->FindNearest(yn,zn,roady,roadz);
3419 Float_t dist = TMath::Sqrt( (yn-cln->GetY())*(yn-cln->GetY())+(zn-cln->GetZ())*(zn-cln->GetZ()));
3422 erry = (dist+0.3)*cln->GetSigmaY2()/TMath::Sqrt(cln->GetQ())*(1.+1./(ddrow));
3423 errz = (dist+0.3)*cln->GetSigmaZ2()/TMath::Sqrt(cln->GetQ())*(1.+1./(ddrow));
3424 if (cln->IsUsed(10)) {
3425 // printf("used\n");
3433 polytrack.AddPoint(xn,cln->GetY(),cln->GetZ(),erry, errz);
3438 if ( (sumused>3) || (sumused>0.5*nfound) || (nfound<0.6*nfoundable)) break;
3439 polytrack.UpdateParameters();
3442 if ( (sumused>3) || (sumused>0.5*nfound)) {
3443 //printf("sumused %d\n",sumused);
3448 polytrack.GetFitDerivation(kr0.GetX(),dy,dz);
3449 AliTPCpolyTrack track2;
3451 polytrack.Refit(track2,0.5+TMath::Abs(dy)*0.3,0.4+TMath::Abs(dz)*0.3);
3452 if (track2.GetN()<0.5*nfoundable) continue;
3455 if ((nfound>0.6*nfoundable) &&( nfoundable>0.4*(i1-i2))) {
3457 // test seed with and without constrain
3458 for (Int_t constrain=0; constrain<=0;constrain++){
3459 // add polytrack candidate
3461 Double_t x[5], c[15];
3462 Double_t x1,x2,x3,y1,y2,y3,z1,z2,z3;
3463 track2.GetBoundaries(x3,x1);
3465 track2.GetFitPoint(x1,y1,z1);
3466 track2.GetFitPoint(x2,y2,z2);
3467 track2.GetFitPoint(x3,y3,z3);
3469 //is track pointing to the vertex ?
3472 polytrack.GetFitPoint(x0,y0,z0);
3485 x[4]=F1(x1,y1,x2,y2,x3,y3);
3487 // if (TMath::Abs(x[4]) >= cuts[0]) continue; //
3488 x[2]=F2(x1,y1,x2,y2,x3,y3);
3490 //if (TMath::Abs(x[4]*x1-x[2]) >= cuts[1]) continue;
3491 //x[3]=F3(x1,y1,x2,y2,z1,z2);
3492 x[3]=F3n(x1,y1,x3,y3,z1,z3,x[4]);
3493 //if (TMath::Abs(x[3]) > cuts[3]) continue;
3496 Double_t sy =0.1, sz =0.1;
3497 Double_t sy1=0.02, sz1=0.02;
3498 Double_t sy2=0.02, sz2=0.02;
3502 sy3=25000*x[4]*x[4]+0.1, sy=0.1, sz=0.1;
3505 Double_t f40=(F1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
3506 Double_t f42=(F1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
3507 Double_t f43=(F1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
3508 Double_t f20=(F2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
3509 Double_t f22=(F2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
3510 Double_t f23=(F2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
3512 Double_t f30=(F3(x1,y1+sy,x3,y3,z1,z3)-x[3])/sy;
3513 Double_t f31=(F3(x1,y1,x3,y3,z1+sz,z3)-x[3])/sz;
3514 Double_t f32=(F3(x1,y1,x3,y3+sy,z1,z3)-x[3])/sy;
3515 Double_t f34=(F3(x1,y1,x3,y3,z1,z3+sz)-x[3])/sz;
3520 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
3521 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
3522 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
3523 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
3524 c[13]=f30*sy1*f40+f32*sy2*f42;
3525 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
3527 //Int_t row1 = fSectors->GetRowNumber(x1);
3528 Int_t row1 = GetRowNumber(x1);
3532 AliTPCseed *track=new (seed) AliTPCseed(index, x, c, x1, sec*alpha+shift);
3533 track->fIsSeeding = kTRUE;
3534 Int_t rc=FollowProlongation(*track, i2);
3535 if (constrain) track->fBConstrain =1;
3537 track->fBConstrain =0;
3538 track->fLastPoint = row1+fInnerSec->GetNRows(); // first cluster in track position
3539 track->fFirstPoint = track->fLastPoint;
3541 if (rc==0 || track->GetNumberOfClusters()<(i1-i2)*0.5 ||
3542 track->GetNumberOfClusters() < track->fNFoundable*0.6 ||
3543 track->fNShared>0.4*track->GetNumberOfClusters()) {
3546 seed->~AliTPCseed();
3549 arr->AddLast(track);
3550 seed = new AliTPCseed;
3554 } // if accepted seed
3557 Info("MakeSeeds2","\nSeeding statiistic:\t%d\t%d\t%d\t%d",nin0,nin1,nin2,nin3);
3563 AliTPCseed *AliTPCtrackerMI::MakeSeed(AliTPCseed *track, Float_t r0, Float_t r1, Float_t r2)
3567 //reseed using track points
3568 Int_t p0 = int(r0*track->GetNumberOfClusters()); // point 0
3569 Int_t p1 = int(r1*track->GetNumberOfClusters());
3570 Int_t p2 = int(r2*track->GetNumberOfClusters()); // last point
3572 Double_t x0[3],x1[3],x2[3];
3573 for (Int_t i=0;i<3;i++){
3579 // find track position at given ratio of the length
3580 Int_t sec0=0, sec1=0, sec2=0;
3583 for (Int_t i=0;i<160;i++){
3584 if (track->fClusterPointer[i]){
3586 AliTPCTrackerPoint *trpoint =track->GetTrackPoint(i);
3587 if ( (index<p0) || x0[0]<0 ){
3588 if (trpoint->GetX()>1){
3589 clindex = track->GetClusterIndex2(i);
3591 x0[0] = trpoint->GetX();
3592 x0[1] = trpoint->GetY();
3593 x0[2] = trpoint->GetZ();
3594 sec0 = ((clindex&0xff000000)>>24)%18;
3599 if ( (index<p1) &&(trpoint->GetX()>1)){
3600 clindex = track->GetClusterIndex2(i);
3602 x1[0] = trpoint->GetX();
3603 x1[1] = trpoint->GetY();
3604 x1[2] = trpoint->GetZ();
3605 sec1 = ((clindex&0xff000000)>>24)%18;
3608 if ( (index<p2) &&(trpoint->GetX()>1)){
3609 clindex = track->GetClusterIndex2(i);
3611 x2[0] = trpoint->GetX();
3612 x2[1] = trpoint->GetY();
3613 x2[2] = trpoint->GetZ();
3614 sec2 = ((clindex&0xff000000)>>24)%18;
3621 Double_t alpha, cs,sn, xx2,yy2;