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 //-------------------------------------------------------
35 #include "Riostream.h"
36 #include <TClonesArray.h>
38 #include <TObjArray.h>
41 #include "AliComplexCluster.h"
44 #include "AliRunLoader.h"
45 #include "AliTPCClustersRow.h"
46 #include "AliTPCParam.h"
47 #include "AliTPCclusterMI.h"
48 #include "AliTPCpolyTrack.h"
49 #include "AliTPCreco.h"
50 #include "AliTPCtrackerMI.h"
51 #include "TStopwatch.h"
53 #include "AliTPCReconstructor.h"
57 ClassImp(AliTPCtrackerMI)
60 class AliTPCFastMath {
63 static Double_t FastAsin(Double_t x);
65 static Double_t fgFastAsin[20000]; //lookup table for fast asin computation
68 Double_t AliTPCFastMath::fgFastAsin[20000];
69 AliTPCFastMath gAliTPCFastMath; // needed to fill the LUT
71 AliTPCFastMath::AliTPCFastMath(){
73 // initialized lookup table;
74 for (Int_t i=0;i<10000;i++){
75 fgFastAsin[2*i] = TMath::ASin(i/10000.);
76 fgFastAsin[2*i+1] = (TMath::ASin((i+1)/10000.)-fgFastAsin[2*i]);
80 Double_t AliTPCFastMath::FastAsin(Double_t x){
82 // return asin using lookup table
84 Int_t index = int(x*10000);
85 return fgFastAsin[2*index]+(x*10000.-index)*fgFastAsin[2*index+1];
88 Int_t index = int(x*10000);
89 return -(fgFastAsin[2*index]+(x*10000.-index)*fgFastAsin[2*index+1]);
95 Int_t AliTPCtrackerMI::UpdateTrack(AliTPCseed * track, Int_t accept){
97 //update track information using current cluster - track->fCurrentCluster
100 AliTPCclusterMI* c =track->fCurrentCluster;
101 if (accept>0) track->fCurrentClusterIndex1 |=0x8000; //sign not accepted clusters
103 UInt_t i = track->fCurrentClusterIndex1;
105 Int_t sec=(i&0xff000000)>>24;
106 //Int_t row = (i&0x00ff0000)>>16;
107 track->fRow=(i&0x00ff0000)>>16;
108 track->fSector = sec;
109 // Int_t index = i&0xFFFF;
110 if (sec>=fParam->GetNInnerSector()) track->fRow += fParam->GetNRowLow();
111 track->SetClusterIndex2(track->fRow, i);
112 //track->fFirstPoint = row;
113 //if ( track->fLastPoint<row) track->fLastPoint =row;
114 // if (track->fRow<0 || track->fRow>160) {
115 // printf("problem\n");
117 if (track->fFirstPoint>track->fRow)
118 track->fFirstPoint = track->fRow;
119 if (track->fLastPoint<track->fRow)
120 track->fLastPoint = track->fRow;
123 track->fClusterPointer[track->fRow] = c;
126 Float_t angle2 = track->GetSnp()*track->GetSnp();
127 angle2 = TMath::Sqrt(angle2/(1-angle2));
129 //SET NEW Track Point
133 AliTPCTrackerPoint &point =*(track->GetTrackPoint(track->fRow));
135 point.SetSigmaY(c->GetSigmaY2()/track->fCurrentSigmaY2);
136 point.SetSigmaZ(c->GetSigmaZ2()/track->fCurrentSigmaZ2);
137 point.SetErrY(sqrt(track->fErrorY2));
138 point.SetErrZ(sqrt(track->fErrorZ2));
140 point.SetX(track->GetX());
141 point.SetY(track->GetY());
142 point.SetZ(track->GetZ());
143 point.SetAngleY(angle2);
144 point.SetAngleZ(track->GetTgl());
145 if (point.fIsShared){
146 track->fErrorY2 *= 4;
147 track->fErrorZ2 *= 4;
151 Double_t chi2 = track->GetPredictedChi2(track->fCurrentCluster);
153 track->fErrorY2 *= 1.3;
154 track->fErrorY2 += 0.01;
155 track->fErrorZ2 *= 1.3;
156 track->fErrorZ2 += 0.005;
158 if (accept>0) return 0;
159 if (track->GetNumberOfClusters()%20==0){
160 // if (track->fHelixIn){
161 // TClonesArray & larr = *(track->fHelixIn);
162 // Int_t ihelix = larr.GetEntriesFast();
163 // new(larr[ihelix]) AliHelix(*track) ;
166 track->fNoCluster =0;
167 return track->Update(c,chi2,i);
172 Int_t AliTPCtrackerMI::AcceptCluster(AliTPCseed * seed, AliTPCclusterMI * cluster, Float_t factor,
173 Float_t cory, Float_t corz)
176 // decide according desired precision to accept given
177 // cluster for tracking
178 Double_t sy2=ErrY2(seed,cluster)*cory;
179 Double_t sz2=ErrZ2(seed,cluster)*corz;
180 //sy2=ErrY2(seed,cluster)*cory;
181 //sz2=ErrZ2(seed,cluster)*cory;
183 Double_t sdistancey2 = sy2+seed->GetSigmaY2();
184 Double_t sdistancez2 = sz2+seed->GetSigmaZ2();
186 Double_t rdistancey2 = (seed->fCurrentCluster->GetY()-seed->GetY())*
187 (seed->fCurrentCluster->GetY()-seed->GetY())/sdistancey2;
188 Double_t rdistancez2 = (seed->fCurrentCluster->GetZ()-seed->GetZ())*
189 (seed->fCurrentCluster->GetZ()-seed->GetZ())/sdistancez2;
191 Double_t rdistance2 = rdistancey2+rdistancez2;
194 if (rdistance2>16) return 3;
197 if ((rdistancey2>9.*factor || rdistancez2>9.*factor) && cluster->GetType()==0)
198 return 2; //suspisiouce - will be changed
200 if ((rdistancey2>6.25*factor || rdistancez2>6.25*factor) && cluster->GetType()>0)
201 // strict cut on overlaped cluster
202 return 2; //suspisiouce - will be changed
204 if ( (rdistancey2>1.*factor || rdistancez2>6.25*factor )
205 && cluster->GetType()<0){
215 //_____________________________________________________________________________
216 AliTPCtrackerMI::AliTPCtrackerMI(const AliTPCParam *par):
217 AliTracker(), fkNIS(par->GetNInnerSector()/2), fkNOS(par->GetNOuterSector()/2)
219 //---------------------------------------------------------------------
220 // The main TPC tracker constructor
221 //---------------------------------------------------------------------
222 fInnerSec=new AliTPCSector[fkNIS];
223 fOuterSec=new AliTPCSector[fkNOS];
226 for (i=0; i<fkNIS; i++) fInnerSec[i].Setup(par,0);
227 for (i=0; i<fkNOS; i++) fOuterSec[i].Setup(par,1);
234 Int_t nrowlow = par->GetNRowLow();
235 Int_t nrowup = par->GetNRowUp();
238 for (Int_t i=0;i<nrowlow;i++){
239 fXRow[i] = par->GetPadRowRadiiLow(i);
240 fPadLength[i]= par->GetPadPitchLength(0,i);
241 fYMax[i] = fXRow[i]*TMath::Tan(0.5*par->GetInnerAngle());
245 for (Int_t i=0;i<nrowup;i++){
246 fXRow[i+nrowlow] = par->GetPadRowRadiiUp(i);
247 fPadLength[i+nrowlow] = par->GetPadPitchLength(60,i);
248 fYMax[i+nrowlow] = fXRow[i+nrowlow]*TMath::Tan(0.5*par->GetOuterAngle());
260 //________________________________________________________________________
261 AliTPCtrackerMI::AliTPCtrackerMI(const AliTPCtrackerMI &t):
266 //------------------------------------
267 // dummy copy constructor
268 //------------------------------------------------------------------
270 AliTPCtrackerMI & AliTPCtrackerMI::operator=(const AliTPCtrackerMI& /*r*/){
271 //------------------------------
273 //--------------------------------------------------------------
276 //_____________________________________________________________________________
277 AliTPCtrackerMI::~AliTPCtrackerMI() {
278 //------------------------------------------------------------------
279 // TPC tracker destructor
280 //------------------------------------------------------------------
289 void AliTPCtrackerMI::SetIO()
293 fInput = AliRunLoader::GetTreeR("TPC", kFALSE,AliConfig::GetDefaultEventFolderName());
295 fOutput = AliRunLoader::GetTreeT("TPC", kTRUE,AliConfig::GetDefaultEventFolderName());
297 AliTPCtrack *iotrack= new AliTPCtrack;
298 fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
304 void AliTPCtrackerMI::SetIO(TTree * input, TTree * output, AliESD * event)
320 AliTPCtrack *iotrack= new AliTPCtrack;
321 // iotrack->fHelixIn = new TClonesArray("AliHelix");
322 //iotrack->fHelixOut = new TClonesArray("AliHelix");
323 fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
326 if (output && (fDebug&2)){
327 //write the full seed information if specified in debug mode
329 fSeedTree = new TTree("Seeds","Seeds");
330 AliTPCseed * vseed = new AliTPCseed;
332 TClonesArray * arrtr = new TClonesArray("AliTPCTrackPoint",160);
333 arrtr->ExpandCreateFast(160);
334 TClonesArray * arre = new TClonesArray("AliTPCExactPoint",160);
336 vseed->fPoints = arrtr;
337 vseed->fEPoints = arre;
338 // vseed->fClusterPoints = arrcl;
339 fSeedTree->Branch("seeds","AliTPCseed",&vseed,32000,99);
342 fTreeDebug = new TTree("trackDebug","trackDebug");
343 TClonesArray * arrd = new TClonesArray("AliTPCTrackPoint2",0);
344 fTreeDebug->Branch("debug",&arrd,32000,99);
352 void AliTPCtrackerMI::FillESD(TObjArray* arr)
356 //fill esds using updated tracks
358 // write tracks to the event
359 // store index of the track
360 Int_t nseed=arr->GetEntriesFast();
361 for (Int_t i=0; i<nseed; i++) {
362 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
364 pt->PropagateTo(fParam->GetInnerRadiusLow());
365 if ( (pt->GetNumberOfClusters()>70)&& (Float_t(pt->GetNumberOfClusters())/Float_t(pt->fNFoundable))>0.55) {
367 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
368 //iotrack.SetTPCindex(i);
369 fEvent->AddTrack(&iotrack);
372 if ( (pt->GetNumberOfClusters()>30) && (Float_t(pt->GetNumberOfClusters())/Float_t(pt->fNFoundable))>0.70) {
373 Int_t found,foundable,shared;
374 pt->GetClusterStatistic(0,60,found, foundable,shared,kFALSE);
375 if ( (found>20) && (pt->fNShared/float(pt->GetNumberOfClusters())<0.2)){
377 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
378 //iotrack.SetTPCindex(i);
379 fEvent->AddTrack(&iotrack);
384 if ( (pt->GetNumberOfClusters()>20) && (Float_t(pt->GetNumberOfClusters())/Float_t(pt->fNFoundable))>0.8) {
385 Int_t found,foundable,shared;
386 pt->GetClusterStatistic(0,60,found, foundable,shared,kFALSE);
387 if (found<20) continue;
388 if (pt->fNShared/float(pt->GetNumberOfClusters())>0.2) continue;
391 iotrack.UpdateTrackParams(pt,AliESDtrack::kTPCin);
392 //iotrack.SetTPCindex(i);
393 fEvent->AddTrack(&iotrack);
400 void AliTPCtrackerMI::WriteTracks(TTree * tree)
403 // write tracks from seed array to selected tree
407 AliTPCtrack *iotrack= new AliTPCtrack;
408 fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
413 void AliTPCtrackerMI::WriteTracks()
416 // write tracks to the given output tree -
417 // output specified with SetIO routine
424 AliTPCtrack *iotrack= 0;
425 Int_t nseed=fSeeds->GetEntriesFast();
426 //for (Int_t i=0; i<nseed; i++) {
427 // iotrack= (AliTPCtrack*)fSeeds->UncheckedAt(i);
428 // if (iotrack) break;
430 //TBranch * br = fOutput->Branch("tracks","AliTPCtrack",&iotrack,32000,100);
431 TBranch * br = fOutput->GetBranch("tracks");
432 br->SetAddress(&iotrack);
434 for (Int_t i=0; i<nseed; i++) {
435 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i);
437 AliTPCtrack * track = new AliTPCtrack(*pt);
440 // br->SetAddress(&iotrack);
445 //fOutput->GetDirectory()->cd();
451 //write the full seed information if specified in debug mode
453 AliTPCseed * vseed = new AliTPCseed;
455 TClonesArray * arrtr = new TClonesArray("AliTPCTrackPoint",160);
456 arrtr->ExpandCreateFast(160);
457 //TClonesArray * arrcl = new TClonesArray("AliTPCclusterMI",160);
458 //arrcl->ExpandCreateFast(160);
459 TClonesArray * arre = new TClonesArray("AliTPCExactPoint",160);
461 vseed->fPoints = arrtr;
462 vseed->fEPoints = arre;
463 // vseed->fClusterPoints = arrcl;
464 //TBranch * brseed = seedtree->Branch("seeds","AliTPCseed",&vseed,32000,99);
465 TBranch * brseed = fSeedTree->GetBranch("seeds");
467 Int_t nseed=fSeeds->GetEntriesFast();
469 for (Int_t i=0; i<nseed; i++) {
470 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i);
473 // pt->fClusterPoints = arrcl;
477 brseed->SetAddress(&vseed);
481 // pt->fClusterPoints = 0;
484 if (fTreeDebug) fTreeDebug->Write();
492 Double_t AliTPCtrackerMI::ErrY2(AliTPCseed* seed, AliTPCclusterMI * cl){
495 //seed->SetErrorY2(0.1);
497 //calculate look-up table at the beginning
498 static Bool_t ginit = kFALSE;
499 static Float_t gnoise1,gnoise2,gnoise3;
500 static Float_t ggg1[10000];
501 static Float_t ggg2[10000];
502 static Float_t ggg3[10000];
503 static Float_t glandau1[10000];
504 static Float_t glandau2[10000];
505 static Float_t glandau3[10000];
507 static Float_t gcor01[500];
508 static Float_t gcor02[500];
509 static Float_t gcorp[500];
514 for (Int_t i=1;i<500;i++){
515 Float_t rsigma = float(i)/100.;
516 gcor02[i] = TMath::Max(0.78 +TMath::Exp(7.4*(rsigma-1.2)),0.6);
517 gcor01[i] = TMath::Max(0.72 +TMath::Exp(3.36*(rsigma-1.2)),0.6);
518 gcorp[i] = TMath::Max(TMath::Power((rsigma+0.5),1.5),1.2);
522 for (Int_t i=3;i<10000;i++){
526 Float_t amp = float(i);
527 Float_t padlength =0.75;
528 gnoise1 = 0.0004/padlength;
529 Float_t nel = 0.268*amp;
530 Float_t nprim = 0.155*amp;
531 ggg1[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.001*nel/(padlength*padlength))/nel;
532 glandau1[i] = (2.+0.12*nprim)*0.5* (2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
533 if (glandau1[i]>1) glandau1[i]=1;
534 glandau1[i]*=padlength*padlength/12.;
538 gnoise2 = 0.0004/padlength;
541 ggg2[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
542 glandau2[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
543 if (glandau2[i]>1) glandau2[i]=1;
544 glandau2[i]*=padlength*padlength/12.;
549 gnoise3 = 0.0004/padlength;
552 ggg3[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
553 glandau3[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
554 if (glandau3[i]>1) glandau3[i]=1;
555 glandau3[i]*=padlength*padlength/12.;
563 Int_t amp = int(TMath::Abs(cl->GetQ()));
565 seed->SetErrorY2(1.);
569 Float_t z = TMath::Abs(fParam->GetZLength()-TMath::Abs(seed->GetZ()));
570 Int_t ctype = cl->GetType();
571 Float_t padlength= GetPadPitchLength(seed->fRow);
572 Float_t angle2 = seed->GetSnp()*seed->GetSnp();
573 angle2 = angle2/(1-angle2);
576 Int_t rsigmay = int(100.*cl->GetSigmaY2()/(seed->fCurrentSigmaY2));
579 if (fSectors==fInnerSec){
581 res = ggg1[amp]*z+glandau1[amp]*angle2;
582 if (ctype==0) res *= gcor01[rsigmay];
585 res*= gcorp[rsigmay];
591 res = ggg2[amp]*z+glandau2[amp]*angle2;
592 if (ctype==0) res *= gcor02[rsigmay];
595 res*= gcorp[rsigmay];
600 res = ggg3[amp]*z+glandau3[amp]*angle2;
601 if (ctype==0) res *= gcor02[rsigmay];
604 res*= gcorp[rsigmay];
611 res*=2.4; // overestimate error 2 times
618 seed->SetErrorY2(res);
626 Double_t AliTPCtrackerMI::ErrZ2(AliTPCseed* seed, AliTPCclusterMI * cl){
629 //seed->SetErrorY2(0.1);
631 //calculate look-up table at the beginning
632 static Bool_t ginit = kFALSE;
633 static Float_t gnoise1,gnoise2,gnoise3;
634 static Float_t ggg1[10000];
635 static Float_t ggg2[10000];
636 static Float_t ggg3[10000];
637 static Float_t glandau1[10000];
638 static Float_t glandau2[10000];
639 static Float_t glandau3[10000];
641 static Float_t gcor01[1000];
642 static Float_t gcor02[1000];
643 static Float_t gcorp[1000];
648 for (Int_t i=1;i<1000;i++){
649 Float_t rsigma = float(i)/100.;
650 gcor02[i] = TMath::Max(0.81 +TMath::Exp(6.8*(rsigma-1.2)),0.6);
651 gcor01[i] = TMath::Max(0.72 +TMath::Exp(2.04*(rsigma-1.2)),0.6);
652 gcorp[i] = TMath::Max(TMath::Power((rsigma+0.5),1.5),1.2);
656 for (Int_t i=3;i<10000;i++){
660 Float_t amp = float(i);
661 Float_t padlength =0.75;
662 gnoise1 = 0.0004/padlength;
663 Float_t nel = 0.268*amp;
664 Float_t nprim = 0.155*amp;
665 ggg1[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.001*nel/(padlength*padlength))/nel;
666 glandau1[i] = (2.+0.12*nprim)*0.5* (2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
667 if (glandau1[i]>1) glandau1[i]=1;
668 glandau1[i]*=padlength*padlength/12.;
672 gnoise2 = 0.0004/padlength;
675 ggg2[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
676 glandau2[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
677 if (glandau2[i]>1) glandau2[i]=1;
678 glandau2[i]*=padlength*padlength/12.;
683 gnoise3 = 0.0004/padlength;
686 ggg3[i] = fParam->GetDiffT()*fParam->GetDiffT()*(2+0.0008*nel/(padlength*padlength))/nel;
687 glandau3[i] = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
688 if (glandau3[i]>1) glandau3[i]=1;
689 glandau3[i]*=padlength*padlength/12.;
697 Int_t amp = int(TMath::Abs(cl->GetQ()));
699 seed->SetErrorY2(1.);
703 Float_t z = TMath::Abs(fParam->GetZLength()-TMath::Abs(seed->GetZ()));
704 Int_t ctype = cl->GetType();
705 Float_t padlength= GetPadPitchLength(seed->fRow);
707 Float_t angle2 = seed->GetSnp()*seed->GetSnp();
708 // if (angle2<0.6) angle2 = 0.6;
709 angle2 = seed->GetTgl()*seed->GetTgl()*(1+angle2/(1-angle2));
712 Int_t rsigmaz = int(100.*cl->GetSigmaZ2()/(seed->fCurrentSigmaZ2));
715 if (fSectors==fInnerSec){
717 res = ggg1[amp]*z+glandau1[amp]*angle2;
718 if (ctype==0) res *= gcor01[rsigmaz];
721 res*= gcorp[rsigmaz];
727 res = ggg2[amp]*z+glandau2[amp]*angle2;
728 if (ctype==0) res *= gcor02[rsigmaz];
731 res*= gcorp[rsigmaz];
736 res = ggg3[amp]*z+glandau3[amp]*angle2;
737 if (ctype==0) res *= gcor02[rsigmaz];
740 res*= gcorp[rsigmaz];
749 if ((ctype<0) &&<70){
757 seed->SetErrorZ2(res);
764 Double_t AliTPCtrackerMI::ErrZ2(AliTPCseed* seed, AliTPCclusterMI * cl){
767 //seed->SetErrorZ2(0.1);
771 Float_t z = TMath::Abs(fParam->GetZLength()-TMath::Abs(seed->GetZ()));
773 Float_t rsigmaz = cl->GetSigmaZ2()/(seed->fCurrentSigmaZ2);
774 Int_t ctype = cl->GetType();
775 Float_t amp = TMath::Abs(cl->GetQ());
780 Float_t landau=2 ; //landau fluctuation part
781 Float_t gg=2; // gg fluctuation part
782 Float_t padlength= GetPadPitchLength(seed->GetX());
784 if (fSectors==fInnerSec){
785 snoise2 = 0.0004/padlength;
788 gg = (2+0.001*nel/(padlength*padlength))/nel;
789 landau = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
790 if (landau>1) landau=1;
793 snoise2 = 0.0004/padlength;
796 gg = (2+0.0008*nel/(padlength*padlength))/nel;
797 landau = (2.+0.12*nprim)*0.5*(2.+nprim*nprim*0.001/(padlength*padlength))/nprim;
798 if (landau>1) landau=1;
800 Float_t sdiff = gg*fParam->GetDiffT()*fParam->GetDiffT()*z;
803 Float_t angle2 = seed->GetSnp()*seed->GetSnp();
804 angle2 = TMath::Sqrt((1-angle2));
805 if (angle2<0.6) angle2 = 0.6;
808 Float_t angle = seed->GetTgl()/angle2;
809 Float_t angular = landau*angle*angle*padlength*padlength/12.;
810 Float_t res = sdiff + angular;
813 if ((ctype==0) && (fSectors ==fOuterSec))
814 res *= 0.81 +TMath::Exp(6.8*(rsigmaz-1.2));
816 if ((ctype==0) && (fSectors ==fInnerSec))
817 res *= 0.72 +TMath::Exp(2.04*(rsigmaz-1.2));
821 res*= TMath::Power(rsigmaz+0.5,1.5); //0.31+0.147*ctype;
827 if ((ctype<0) &&<70){
835 seed->SetErrorZ2(res);
842 void AliTPCseed::Reset(Bool_t all)
846 SetNumberOfClusters(0);
852 for (Int_t i=0;i<8;i++){
853 delete [] fTrackPoints[i];
861 for (Int_t i=0;i<200;i++) SetClusterIndex2(i,-3);
862 for (Int_t i=0;i<160;i++) fClusterPointer[i]=0;
868 void AliTPCseed::Modify(Double_t factor)
871 //------------------------------------------------------------------
872 //This function makes a track forget its history :)
873 //------------------------------------------------------------------
879 fC10*=0; fC11*=factor;
880 fC20*=0; fC21*=0; fC22*=factor;
881 fC30*=0; fC31*=0; fC32*=0; fC33*=factor;
882 fC40*=0; fC41*=0; fC42*=0; fC43*=0; fC44*=factor;
883 SetNumberOfClusters(0);
887 fCurrentSigmaY2 = 0.000005;
888 fCurrentSigmaZ2 = 0.000005;
897 Int_t AliTPCseed::GetProlongation(Double_t xk, Double_t &y, Double_t & z) const
899 //-----------------------------------------------------------------
900 // This function find proloncation of a track to a reference plane x=xk.
901 // doesn't change internal state of the track
902 //-----------------------------------------------------------------
904 Double_t x1=fX, x2=x1+(xk-x1), dx=x2-x1;
906 if (TMath::Abs(fP4*xk - fP2) >= 0.999) {
910 // Double_t y1=fP0, z1=fP1;
911 Double_t c1=fP4*x1 - fP2, r1=sqrt(1.- c1*c1);
912 Double_t c2=fP4*x2 - fP2, r2=sqrt(1.- c2*c2);
916 //y += dx*(c1+c2)/(r1+r2);
917 //z += dx*(c1+c2)/(c1*r2 + c2*r1)*fP3;
919 Double_t dy = dx*(c1+c2)/(r1+r2);
922 Double_t delta = fP4*dx*(c1+c2)/(c1*r2 + c2*r1);
924 if (TMath::Abs(delta)>0.0001){
925 dz = fP3*TMath::ASin(delta)/fP4;
927 dz = dx*fP3*(c1+c2)/(c1*r2 + c2*r1);
930 dz = fP3*AliTPCFastMath::FastAsin(delta)/fP4;
940 //_____________________________________________________________________________
941 Double_t AliTPCseed::GetPredictedChi2(const AliTPCclusterMI *c) const
943 //-----------------------------------------------------------------
944 // This function calculates a predicted chi2 increment.
945 //-----------------------------------------------------------------
946 //Double_t r00=c->GetSigmaY2(), r01=0., r11=c->GetSigmaZ2();
947 Double_t r00=fErrorY2, r01=0., r11=fErrorZ2;
948 r00+=fC00; r01+=fC10; r11+=fC11;
950 Double_t det=r00*r11 - r01*r01;
951 if (TMath::Abs(det) < 1.e-10) {
952 Int_t n=GetNumberOfClusters();
953 if (n>4) cerr<<n<<" AliKalmanTrack warning: Singular matrix !\n";
956 Double_t tmp=r00; r00=r11; r11=tmp; r01=-r01;
958 Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
960 return (dy*r00*dy + 2*r01*dy*dz + dz*r11*dz)/det;
964 //_________________________________________________________________________________________
967 Int_t AliTPCseed::Compare(const TObject *o) const {
968 //-----------------------------------------------------------------
969 // This function compares tracks according to the sector - for given sector according z
970 //-----------------------------------------------------------------
971 AliTPCseed *t=(AliTPCseed*)o;
974 if (t->fRelativeSector>fRelativeSector) return -1;
975 if (t->fRelativeSector<fRelativeSector) return 1;
976 Double_t z2 = t->GetZ();
977 Double_t z1 = GetZ();
979 if (z2<z1) return -1;
984 f2 = 1-20*TMath::Sqrt(t->fC44)/(TMath::Abs(t->GetC())+0.0066);
985 if (t->fBConstrain) f2=1.2;
988 f1 = 1-20*TMath::Sqrt(fC44)/(TMath::Abs(GetC())+0.0066);
990 if (fBConstrain) f1=1.2;
992 if (t->GetNumberOfClusters()*f2 <GetNumberOfClusters()*f1) return -1;
997 void AliTPCtrackerMI::RotateToLocal(AliTPCseed *seed)
999 //rotate to track "local coordinata
1000 Float_t x = seed->GetX();
1001 Float_t y = seed->GetY();
1002 Float_t ymax = x*TMath::Tan(0.5*fSectors->GetAlpha());
1005 seed->fRelativeSector= (seed->fRelativeSector+1) % fN;
1006 if (!seed->Rotate(fSectors->GetAlpha()))
1008 } else if (y <-ymax) {
1009 seed->fRelativeSector= (seed->fRelativeSector-1+fN) % fN;
1010 if (!seed->Rotate(-fSectors->GetAlpha()))
1019 //_____________________________________________________________________________
1020 Int_t AliTPCseed::Update(const AliTPCclusterMI *c, Double_t chisq, UInt_t /*index*/) {
1021 //-----------------------------------------------------------------
1022 // This function associates a cluster with this track.
1023 //-----------------------------------------------------------------
1024 Double_t r00=fErrorY2, r01=0., r11=fErrorZ2;
1026 r00+=fC00; r01+=fC10; r11+=fC11;
1027 Double_t det=r00*r11 - r01*r01;
1028 Double_t tmp=r00; r00=r11/det; r11=tmp/det; r01=-r01/det;
1030 Double_t k00=fC00*r00+fC10*r01, k01=fC00*r01+fC10*r11;
1031 Double_t k10=fC10*r00+fC11*r01, k11=fC10*r01+fC11*r11;
1032 Double_t k20=fC20*r00+fC21*r01, k21=fC20*r01+fC21*r11;
1033 Double_t k30=fC30*r00+fC31*r01, k31=fC30*r01+fC31*r11;
1034 Double_t k40=fC40*r00+fC41*r01, k41=fC40*r01+fC41*r11;
1036 Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
1037 Double_t cur=fP4 + k40*dy + k41*dz, eta=fP2 + k20*dy + k21*dz;
1038 if (TMath::Abs(cur*fX-eta) >= 0.9) {
1042 fP0 += k00*dy + k01*dz;
1043 fP1 += k10*dy + k11*dz;
1045 fP3 += k30*dy + k31*dz;
1048 Double_t c01=fC10, c02=fC20, c03=fC30, c04=fC40;
1049 Double_t c12=fC21, c13=fC31, c14=fC41;
1051 fC00-=k00*fC00+k01*fC10; fC10-=k00*c01+k01*fC11;
1052 fC20-=k00*c02+k01*c12; fC30-=k00*c03+k01*c13;
1053 fC40-=k00*c04+k01*c14;
1055 fC11-=k10*c01+k11*fC11;
1056 fC21-=k10*c02+k11*c12; fC31-=k10*c03+k11*c13;
1057 fC41-=k10*c04+k11*c14;
1059 fC22-=k20*c02+k21*c12; fC32-=k20*c03+k21*c13;
1060 fC42-=k20*c04+k21*c14;
1062 fC33-=k30*c03+k31*c13;
1063 fC43-=k40*c03+k41*c13;
1065 fC44-=k40*c04+k41*c14;
1067 Int_t n=GetNumberOfClusters();
1069 SetNumberOfClusters(n+1);
1070 SetChi2(GetChi2()+chisq);
1077 //_____________________________________________________________________________
1078 Double_t AliTPCtrackerMI::F1old(Double_t x1,Double_t y1,
1079 Double_t x2,Double_t y2,
1080 Double_t x3,Double_t y3)
1082 //-----------------------------------------------------------------
1083 // Initial approximation of the track curvature
1084 //-----------------------------------------------------------------
1085 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
1086 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
1087 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
1088 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
1089 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
1091 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
1092 if ( xr*xr+yr*yr<=0.00000000000001) return 100;
1093 return -xr*yr/sqrt(xr*xr+yr*yr);
1098 //_____________________________________________________________________________
1099 Double_t AliTPCtrackerMI::F1(Double_t x1,Double_t y1,
1100 Double_t x2,Double_t y2,
1101 Double_t x3,Double_t y3)
1103 //-----------------------------------------------------------------
1104 // Initial approximation of the track curvature
1105 //-----------------------------------------------------------------
1111 Double_t det = x3*y2-x2*y3;
1116 Double_t u = 0.5* (x2*(x2-x3)+y2*(y2-y3))/det;
1117 Double_t x0 = x3*0.5-y3*u;
1118 Double_t y0 = y3*0.5+x3*u;
1119 Double_t c2 = 1/TMath::Sqrt(x0*x0+y0*y0);
1125 Double_t AliTPCtrackerMI::F2(Double_t x1,Double_t y1,
1126 Double_t x2,Double_t y2,
1127 Double_t x3,Double_t y3)
1129 //-----------------------------------------------------------------
1130 // Initial approximation of the track curvature
1131 //-----------------------------------------------------------------
1137 Double_t det = x3*y2-x2*y3;
1142 Double_t u = 0.5* (x2*(x2-x3)+y2*(y2-y3))/det;
1143 Double_t x0 = x3*0.5-y3*u;
1144 Double_t y0 = y3*0.5+x3*u;
1145 Double_t c2 = 1/TMath::Sqrt(x0*x0+y0*y0);
1154 //_____________________________________________________________________________
1155 Double_t AliTPCtrackerMI::F2old(Double_t x1,Double_t y1,
1156 Double_t x2,Double_t y2,
1157 Double_t x3,Double_t y3)
1159 //-----------------------------------------------------------------
1160 // Initial approximation of the track curvature times center of curvature
1161 //-----------------------------------------------------------------
1162 Double_t d=(x2-x1)*(y3-y2)-(x3-x2)*(y2-y1);
1163 Double_t a=0.5*((y3-y2)*(y2*y2-y1*y1+x2*x2-x1*x1)-
1164 (y2-y1)*(y3*y3-y2*y2+x3*x3-x2*x2));
1165 Double_t b=0.5*((x2-x1)*(y3*y3-y2*y2+x3*x3-x2*x2)-
1166 (x3-x2)*(y2*y2-y1*y1+x2*x2-x1*x1));
1168 Double_t xr=TMath::Abs(d/(d*x1-a)), yr=d/(d*y1-b);
1170 return -a/(d*y1-b)*xr/sqrt(xr*xr+yr*yr);
1173 //_____________________________________________________________________________
1174 Double_t AliTPCtrackerMI::F3(Double_t x1,Double_t y1,
1175 Double_t x2,Double_t y2,
1176 Double_t z1,Double_t z2)
1178 //-----------------------------------------------------------------
1179 // Initial approximation of the tangent of the track dip angle
1180 //-----------------------------------------------------------------
1181 return (z1 - z2)/sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
1185 Double_t AliTPCtrackerMI::F3n(Double_t x1,Double_t y1,
1186 Double_t x2,Double_t y2,
1187 Double_t z1,Double_t z2, Double_t c)
1189 //-----------------------------------------------------------------
1190 // Initial approximation of the tangent of the track dip angle
1191 //-----------------------------------------------------------------
1195 //angle1 = (z1-z2)*c/(TMath::ASin(c*x1-ni)-TMath::ASin(c*x2-ni));
1197 Double_t d = TMath::Sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2));
1198 if (TMath::Abs(d*c*0.5)>1) return 0;
1199 // Double_t angle2 = TMath::ASin(d*c*0.5);
1200 // Double_t angle2 = AliTPCFastMath::FastAsin(d*c*0.5);
1201 Double_t angle2 = (d*c*0.5>0.1)? TMath::ASin(d*c*0.5): AliTPCFastMath::FastAsin(d*c*0.5);
1203 angle2 = (z1-z2)*c/(angle2*2.);
1207 Bool_t AliTPCtrackerMI::GetProlongation(Double_t x1, Double_t x2, Double_t x[5], Double_t &y, Double_t &z)
1208 {//-----------------------------------------------------------------
1209 // This function find proloncation of a track to a reference plane x=x2.
1210 //-----------------------------------------------------------------
1214 if (TMath::Abs(x[4]*x1 - x[2]) >= 0.999) {
1218 Double_t c1=x[4]*x1 - x[2], r1=sqrt(1.- c1*c1);
1219 Double_t c2=x[4]*x2 - x[2], r2=sqrt(1.- c2*c2);
1223 Double_t dy = dx*(c1+c2)/(r1+r2);
1226 Double_t delta = x[4]*dx*(c1+c2)/(c1*r2 + c2*r1);
1228 if (TMath::Abs(delta)>0.01){
1229 dz = x[3]*TMath::ASin(delta)/x[4];
1231 dz = x[3]*AliTPCFastMath::FastAsin(delta)/x[4];
1234 //dz = x[3]*AliTPCFastMath::FastAsin(delta)/x[4];
1242 Int_t AliTPCtrackerMI::LoadClusters (TTree *tree)
1247 return LoadClusters();
1250 Int_t AliTPCtrackerMI::LoadClusters()
1253 // load clusters to the memory
1254 AliTPCClustersRow *clrow= new AliTPCClustersRow;
1255 clrow->SetClass("AliTPCclusterMI");
1257 clrow->GetArray()->ExpandCreateFast(10000);
1259 // TTree * tree = fClustersArray.GetTree();
1261 TTree * tree = fInput;
1262 TBranch * br = tree->GetBranch("Segment");
1263 br->SetAddress(&clrow);
1265 Int_t j=Int_t(tree->GetEntries());
1266 for (Int_t i=0; i<j; i++) {
1270 fParam->AdjustSectorRow(clrow->GetID(),sec,row);
1272 AliTPCRow * tpcrow=0;
1275 tpcrow = &(fInnerSec[sec%fkNIS][row]);
1279 tpcrow = &(fOuterSec[(sec-fkNIS*2)%fkNOS][row]);
1280 left = (sec-fkNIS*2)/fkNOS;
1283 tpcrow->fN1 = clrow->GetArray()->GetEntriesFast();
1284 tpcrow->fClusters1 = new AliTPCclusterMI[tpcrow->fN1];
1285 for (Int_t i=0;i<tpcrow->fN1;i++)
1286 tpcrow->fClusters1[i] = *(AliTPCclusterMI*)(clrow->GetArray()->At(i));
1289 tpcrow->fN2 = clrow->GetArray()->GetEntriesFast();
1290 tpcrow->fClusters2 = new AliTPCclusterMI[tpcrow->fN2];
1291 for (Int_t i=0;i<tpcrow->fN2;i++)
1292 tpcrow->fClusters2[i] = *(AliTPCclusterMI*)(clrow->GetArray()->At(i));
1303 void AliTPCtrackerMI::UnloadClusters()
1306 // unload clusters from the memory
1308 Int_t nrows = fOuterSec->GetNRows();
1309 for (Int_t sec = 0;sec<fkNOS;sec++)
1310 for (Int_t row = 0;row<nrows;row++){
1311 AliTPCRow* tpcrow = &(fOuterSec[sec%fkNOS][row]);
1313 // if (tpcrow->fClusters1) delete []tpcrow->fClusters1;
1314 // if (tpcrow->fClusters2) delete []tpcrow->fClusters2;
1316 tpcrow->ResetClusters();
1319 nrows = fInnerSec->GetNRows();
1320 for (Int_t sec = 0;sec<fkNIS;sec++)
1321 for (Int_t row = 0;row<nrows;row++){
1322 AliTPCRow* tpcrow = &(fInnerSec[sec%fkNIS][row]);
1324 // if (tpcrow->fClusters1) delete []tpcrow->fClusters1;
1325 //if (tpcrow->fClusters2) delete []tpcrow->fClusters2;
1327 tpcrow->ResetClusters();
1334 //_____________________________________________________________________________
1335 Int_t AliTPCtrackerMI::LoadOuterSectors() {
1336 //-----------------------------------------------------------------
1337 // This function fills outer TPC sectors with clusters.
1338 //-----------------------------------------------------------------
1339 Int_t nrows = fOuterSec->GetNRows();
1341 for (Int_t sec = 0;sec<fkNOS;sec++)
1342 for (Int_t row = 0;row<nrows;row++){
1343 AliTPCRow* tpcrow = &(fOuterSec[sec%fkNOS][row]);
1344 Int_t sec2 = sec+2*fkNIS;
1346 Int_t ncl = tpcrow->fN1;
1348 AliTPCclusterMI *c= &(tpcrow->fClusters1[ncl]);
1349 index=(((sec2<<8)+row)<<16)+ncl;
1350 tpcrow->InsertCluster(c,index);
1355 AliTPCclusterMI *c= &(tpcrow->fClusters2[ncl]);
1356 index=((((sec2+fkNOS)<<8)+row)<<16)+ncl;
1357 tpcrow->InsertCluster(c,index);
1360 // write indexes for fast acces
1362 for (Int_t i=0;i<510;i++)
1363 tpcrow->fFastCluster[i]=-1;
1364 for (Int_t i=0;i<tpcrow->GetN();i++){
1365 Int_t zi = Int_t((*tpcrow)[i]->GetZ()+255.);
1366 tpcrow->fFastCluster[zi]=i; // write index
1369 for (Int_t i=0;i<510;i++){
1370 if (tpcrow->fFastCluster[i]<0)
1371 tpcrow->fFastCluster[i] = last;
1373 last = tpcrow->fFastCluster[i];
1382 //_____________________________________________________________________________
1383 Int_t AliTPCtrackerMI::LoadInnerSectors() {
1384 //-----------------------------------------------------------------
1385 // This function fills inner TPC sectors with clusters.
1386 //-----------------------------------------------------------------
1387 Int_t nrows = fInnerSec->GetNRows();
1389 for (Int_t sec = 0;sec<fkNIS;sec++)
1390 for (Int_t row = 0;row<nrows;row++){
1391 AliTPCRow* tpcrow = &(fInnerSec[sec%fkNIS][row]);
1394 Int_t ncl = tpcrow->fN1;
1396 AliTPCclusterMI *c= &(tpcrow->fClusters1[ncl]);
1397 index=(((sec<<8)+row)<<16)+ncl;
1398 tpcrow->InsertCluster(c,index);
1403 AliTPCclusterMI *c= &(tpcrow->fClusters2[ncl]);
1404 index=((((sec+fkNIS)<<8)+row)<<16)+ncl;
1405 tpcrow->InsertCluster(c,index);
1408 // write indexes for fast acces
1410 for (Int_t i=0;i<510;i++)
1411 tpcrow->fFastCluster[i]=-1;
1412 for (Int_t i=0;i<tpcrow->GetN();i++){
1413 Int_t zi = Int_t((*tpcrow)[i]->GetZ()+255.);
1414 tpcrow->fFastCluster[zi]=i; // write index
1417 for (Int_t i=0;i<510;i++){
1418 if (tpcrow->fFastCluster[i]<0)
1419 tpcrow->fFastCluster[i] = last;
1421 last = tpcrow->fFastCluster[i];
1433 //_________________________________________________________________________
1434 AliTPCclusterMI *AliTPCtrackerMI::GetClusterMI(Int_t index) const {
1435 //--------------------------------------------------------------------
1436 // Return pointer to a given cluster
1437 //--------------------------------------------------------------------
1438 Int_t sec=(index&0xff000000)>>24;
1439 Int_t row=(index&0x00ff0000)>>16;
1440 Int_t ncl=(index&0x00007fff)>>00;
1442 const AliTPCRow * tpcrow=0;
1443 AliTPCclusterMI * clrow =0;
1446 tpcrow = &(fInnerSec[sec%fkNIS][row]);
1447 if (tpcrow==0) return 0;
1450 if (tpcrow->fN1<=ncl) return 0;
1451 clrow = tpcrow->fClusters1;
1454 if (tpcrow->fN2<=ncl) return 0;
1455 clrow = tpcrow->fClusters2;
1459 tpcrow = &(fOuterSec[(sec-fkNIS*2)%fkNOS][row]);
1460 if (tpcrow==0) return 0;
1462 if (sec-2*fkNIS<fkNOS) {
1463 if (tpcrow->fN1<=ncl) return 0;
1464 clrow = tpcrow->fClusters1;
1467 if (tpcrow->fN2<=ncl) return 0;
1468 clrow = tpcrow->fClusters2;
1472 return &(clrow[ncl]);
1478 Int_t AliTPCtrackerMI::FollowToNext(AliTPCseed& t, Int_t nr) {
1479 //-----------------------------------------------------------------
1480 // This function tries to find a track prolongation to next pad row
1481 //-----------------------------------------------------------------
1483 Double_t x= GetXrow(nr), ymax=GetMaxY(nr);
1484 AliTPCclusterMI *cl=0;
1485 Int_t tpcindex= t.GetClusterIndex2(nr);
1487 // update current shape info every 5 pad-row
1488 // if ( (nr%5==0) || t.GetNumberOfClusters()<2 || (t.fCurrentSigmaY2<0.0001) ){
1492 if (fIteration>0 && tpcindex>=-1){ //if we have already clusters
1494 if (tpcindex==-1) return 0; //track in dead zone
1496 cl = t.fClusterPointer[nr];
1497 if ( (cl==0) ) cl = GetClusterMI(tpcindex);
1498 t.fCurrentClusterIndex1 = tpcindex;
1501 Int_t relativesector = ((tpcindex&0xff000000)>>24)%18; // if previously accepted cluster in different sector
1502 Float_t angle = relativesector*fSectors->GetAlpha()+fSectors->GetAlphaShift();
1504 if (angle<-TMath::Pi()) angle += 2*TMath::Pi();
1505 if (angle>=TMath::Pi()) angle -= 2*TMath::Pi();
1507 if (TMath::Abs(angle-t.GetAlpha())>0.001){
1508 Double_t rotation = angle-t.GetAlpha();
1509 t.fRelativeSector= relativesector;
1514 t.fCurrentCluster = cl;
1516 Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1517 if ((tpcindex&0x8000)==0) accept =0;
1519 //if founded cluster is acceptible
1520 if (cl->IsUsed(11)) { // id cluster is shared inrease uncertainty
1527 UpdateTrack(&t,accept);
1532 if (fIteration>1) return 0; // not look for new cluster during refitting
1535 if (TMath::Abs(t.GetSnp())>0.95 || TMath::Abs(x*t.GetC()-t.GetEta())>0.95) return 0;
1536 Double_t y=t.GetYat(x);
1537 if (TMath::Abs(y)>ymax){
1539 t.fRelativeSector= (t.fRelativeSector+1) % fN;
1540 if (!t.Rotate(fSectors->GetAlpha()))
1542 } else if (y <-ymax) {
1543 t.fRelativeSector= (t.fRelativeSector-1+fN) % fN;
1544 if (!t.Rotate(-fSectors->GetAlpha()))
1550 if (!t.PropagateTo(x)) {
1551 if (fIteration==0) t.fRemoval = 10;
1555 Double_t z=t.GetZ();
1557 const AliTPCRow &krow=GetRow(t.fRelativeSector,nr);
1558 if ( (t.GetSigmaY2()<0) || t.GetSigmaZ2()<0) return 0;
1560 Double_t roadz = 1.;
1562 if (TMath::Abs(TMath::Abs(y)-ymax)<krow.fDeadZone){
1564 t.SetClusterIndex2(nr,-1);
1569 if (TMath::Abs(z)<(AliTPCReconstructor::GetCtgRange()*x+10)) t.fNFoundable++;
1575 // cl = krow.FindNearest2(y+10.,z,roady,roadz,index);
1576 cl = krow.FindNearest2(y,z,roady,roadz,index);
1577 if (cl) t.fCurrentClusterIndex1 = krow.GetIndex(index);
1580 t.fCurrentCluster = cl;
1582 if (fIteration==2&&cl->IsUsed(10)) return 0;
1583 Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1584 if (fIteration==2&&cl->IsUsed(11)) {
1591 if (t.fCurrentCluster->IsUsed(10)){
1596 if (t.fNShared>0.7*t.GetNumberOfClusters()) {
1602 if (accept<3) UpdateTrack(&t,accept);
1605 if ( fIteration==0 && t.fNFoundable*0.5 > t.GetNumberOfClusters()) t.fRemoval=10;
1611 Int_t AliTPCtrackerMI::FollowToNextFast(AliTPCseed& t, Int_t nr) {
1612 //-----------------------------------------------------------------
1613 // This function tries to find a track prolongation to next pad row
1614 //-----------------------------------------------------------------
1616 Double_t x= GetXrow(nr), ymax=GetMaxY(nr);
1618 if (!t.GetProlongation(x,y,z)) {
1624 if (TMath::Abs(y)>ymax){
1627 t.fRelativeSector= (t.fRelativeSector+1) % fN;
1628 if (!t.Rotate(fSectors->GetAlpha()))
1630 } else if (y <-ymax) {
1631 t.fRelativeSector= (t.fRelativeSector-1+fN) % fN;
1632 if (!t.Rotate(-fSectors->GetAlpha()))
1635 if (!t.PropagateTo(x)) {
1638 t.GetProlongation(x,y,z);
1641 // update current shape info every 3 pad-row
1642 if ( (nr%6==0) || t.GetNumberOfClusters()<2 || (t.fCurrentSigmaY2<0.0001) ){
1643 // t.fCurrentSigmaY = GetSigmaY(&t);
1644 //t.fCurrentSigmaZ = GetSigmaZ(&t);
1648 AliTPCclusterMI *cl=0;
1653 const AliTPCRow &krow=GetRow(t.fRelativeSector,nr);
1654 if ( (t.GetSigmaY2()<0) || t.GetSigmaZ2()<0) return 0;
1656 Double_t roadz = 1.;
1659 if (TMath::Abs(TMath::Abs(y)-ymax)<krow.fDeadZone){
1661 t.SetClusterIndex2(row,-1);
1666 if (TMath::Abs(z)>(AliTPCReconstructor::GetCtgRange()*x+10)) t.SetClusterIndex2(row,-1);
1670 if ((cl==0)&&(krow)) {
1671 // cl = krow.FindNearest2(y+10,z,roady,roadz,index);
1672 cl = krow.FindNearest2(y,z,roady,roadz,index);
1674 if (cl) t.fCurrentClusterIndex1 = krow.GetIndex(index);
1678 t.fCurrentCluster = cl;
1679 // Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1681 t.SetClusterIndex2(row,index);
1682 t.fClusterPointer[row] = cl;
1690 Int_t AliTPCtrackerMI::UpdateClusters(AliTPCseed& t, Int_t nr) {
1691 //-----------------------------------------------------------------
1692 // This function tries to find a track prolongation to next pad row
1693 //-----------------------------------------------------------------
1694 t.fCurrentCluster = 0;
1695 t.fCurrentClusterIndex1 = 0;
1697 Double_t xt=t.GetX();
1698 Int_t row = GetRowNumber(xt)-1;
1699 Double_t ymax= GetMaxY(nr);
1701 if (row < nr) return 1; // don't prolongate if not information until now -
1702 if (TMath::Abs(t.GetSnp())>0.9 && t.GetNumberOfClusters()>40. && fIteration!=2) {
1704 return 0; // not prolongate strongly inclined tracks
1706 if (TMath::Abs(t.GetSnp())>0.95) {
1708 return 0; // not prolongate strongly inclined tracks
1711 Double_t x= GetXrow(nr);
1713 //t.PropagateTo(x+0.02);
1714 //t.PropagateTo(x+0.01);
1715 if (!t.PropagateTo(x)){
1722 if (TMath::Abs(y)>ymax){
1724 t.fRelativeSector= (t.fRelativeSector+1) % fN;
1725 if (!t.Rotate(fSectors->GetAlpha()))
1727 } else if (y <-ymax) {
1728 t.fRelativeSector= (t.fRelativeSector-1+fN) % fN;
1729 if (!t.Rotate(-fSectors->GetAlpha()))
1732 // if (!t.PropagateTo(x)){
1740 AliTPCRow &krow=GetRow(t.fRelativeSector,nr);
1742 if (TMath::Abs(TMath::Abs(y)-ymax)<krow.fDeadZone){
1744 t.SetClusterIndex2(nr,-1);
1749 if (TMath::Abs(t.GetZ())<(AliTPCReconstructor::GetCtgRange()*t.GetX()+10)) t.fNFoundable++;
1755 if ( (nr%6==0) || t.GetNumberOfClusters()<2){
1756 // t.fCurrentSigmaY = GetSigmaY(&t);
1757 //t.fCurrentSigmaZ = GetSigmaZ(&t);
1761 AliTPCclusterMI *cl=0;
1764 Double_t roady = 1.;
1765 Double_t roadz = 1.;
1769 index = t.GetClusterIndex2(nr);
1770 if ( (index>0) && (index&0x8000)==0){
1771 cl = t.fClusterPointer[nr];
1772 if ( (cl==0) && (index>0)) cl = GetClusterMI(index);
1773 t.fCurrentClusterIndex1 = index;
1775 t.fCurrentCluster = cl;
1782 //cl = krow.FindNearest2(y+10,z,roady,roadz,index);
1783 cl = krow.FindNearest2(y,z,roady,roadz,index);
1786 if (cl) t.fCurrentClusterIndex1 = krow.GetIndex(index);
1787 t.fCurrentCluster = cl;
1793 Int_t AliTPCtrackerMI::FollowToNextCluster(AliTPCseed & t, Int_t nr) {
1794 //-----------------------------------------------------------------
1795 // This function tries to find a track prolongation to next pad row
1796 //-----------------------------------------------------------------
1798 //update error according neighborhoud
1800 if (t.fCurrentCluster) {
1802 Int_t accept = AcceptCluster(&t,t.fCurrentCluster,1.);
1804 if (t.fCurrentCluster->IsUsed(10)){
1810 if (t.fNShared>0.7*t.GetNumberOfClusters()) {
1815 if (fIteration>0) accept = 0;
1816 if (accept<3) UpdateTrack(&t,accept);
1820 if ( ( (t.GetSigmaY2()+t.GetSigmaZ2())>0.16)&& t.GetNumberOfClusters()>18) t.fRemoval=10;
1821 if ( t.GetChi2()/t.GetNumberOfClusters()>6 &&t.GetNumberOfClusters()>18) t.fRemoval=10;
1823 if (( (t.fNFoundable*0.5 > t.GetNumberOfClusters()) || t.fNoCluster>15)) t.fRemoval=10;
1831 //_____________________________________________________________________________
1832 Int_t AliTPCtrackerMI::FollowProlongation(AliTPCseed& t, Int_t rf, Int_t step) {
1833 //-----------------------------------------------------------------
1834 // This function tries to find a track prolongation.
1835 //-----------------------------------------------------------------
1836 Double_t xt=t.GetX();
1838 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
1839 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1840 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1842 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
1844 Int_t first = GetRowNumber(xt)-1;
1845 for (Int_t nr= first; nr>=rf; nr-=step) {
1846 if (nr<fInnerSec->GetNRows())
1847 fSectors = fInnerSec;
1849 fSectors = fOuterSec;
1850 if (FollowToNext(t,nr)==0)
1859 //_____________________________________________________________________________
1860 Int_t AliTPCtrackerMI::FollowProlongationFast(AliTPCseed& t, Int_t rf, Int_t step) {
1861 //-----------------------------------------------------------------
1862 // This function tries to find a track prolongation.
1863 //-----------------------------------------------------------------
1864 Double_t xt=t.GetX();
1866 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
1867 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1868 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1869 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
1871 for (Int_t nr=GetRowNumber(xt)-1; nr>=rf; nr-=step) {
1873 if (FollowToNextFast(t,nr)==0)
1874 if (!t.IsActive()) return 0;
1884 Int_t AliTPCtrackerMI::FollowBackProlongation(AliTPCseed& t, Int_t rf) {
1885 //-----------------------------------------------------------------
1886 // This function tries to find a track prolongation.
1887 //-----------------------------------------------------------------
1888 // Double_t xt=t.GetX();
1890 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
1891 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
1892 if (alpha < 0. ) alpha += 2.*TMath::Pi();
1893 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
1895 Int_t first = t.fFirstPoint;
1897 if (first<0) first=0;
1898 for (Int_t nr=first; nr<=rf; nr++) {
1899 //if ( (t.GetSnp()<0.9))
1900 if (nr<fInnerSec->GetNRows())
1901 fSectors = fInnerSec;
1903 fSectors = fOuterSec;
1913 Float_t AliTPCtrackerMI::OverlapFactor(AliTPCseed * s1, AliTPCseed * s2, Int_t &sum1, Int_t & sum2)
1921 Float_t dz2 =(s1->GetZ() - s2->GetZ());
1924 Float_t dy2 =TMath::Abs((s1->GetY() - s2->GetY()));
1926 Float_t distance = TMath::Sqrt(dz2+dy2);
1927 if (distance>4.) return 0; // if there are far away - not overlap - to reduce combinatorics
1930 Int_t firstpoint = TMath::Min(s1->fFirstPoint,s2->fFirstPoint);
1931 Int_t lastpoint = TMath::Max(s1->fLastPoint,s2->fLastPoint);
1936 if (firstpoint>lastpoint) {
1937 firstpoint =lastpoint;
1942 for (Int_t i=firstpoint-1;i<lastpoint+1;i++){
1943 if (s1->GetClusterIndex2(i)>0) sum1++;
1944 if (s2->GetClusterIndex2(i)>0) sum2++;
1945 if (s1->GetClusterIndex2(i)==s2->GetClusterIndex2(i) && s1->GetClusterIndex2(i)>0) {
1949 if (sum<5) return 0;
1951 Float_t summin = TMath::Min(sum1+1,sum2+1);
1952 Float_t ratio = (sum+1)/Float_t(summin);
1956 void AliTPCtrackerMI::SignShared(AliTPCseed * s1, AliTPCseed * s2)
1960 if (TMath::Abs(s1->GetC()-s2->GetC())>0.004) return;
1961 if (TMath::Abs(s1->GetTgl()-s2->GetTgl())>0.6) return;
1963 Float_t dz2 =(s1->GetZ() - s2->GetZ());
1965 Float_t dy2 =(s1->GetY() - s2->GetY());
1967 Float_t distance = dz2+dy2;
1968 if (distance>325.) return ; // if there are far away - not overlap - to reduce combinatorics
1973 Int_t firstpoint = TMath::Max(s1->fFirstPoint,s2->fFirstPoint);
1974 Int_t lastpoint = TMath::Min(s1->fLastPoint,s2->fLastPoint);
1976 if (firstpoint>=lastpoint-5) return;;
1978 for (Int_t i=firstpoint;i<lastpoint;i++){
1979 // if ( (s1->GetClusterIndex2(i)&0xFFFF8FFF)==(s2->GetClusterIndex2(i)&0xFFFF8FFF) && s1->GetClusterIndex2(i)>0) {
1980 if ( (s1->GetClusterIndex2(i))==(s2->GetClusterIndex2(i)) && s1->GetClusterIndex2(i)>0) {
1987 for (Int_t i=firstpoint;i<lastpoint;i++){
1988 // if ( (s1->GetClusterIndex2(i)&0xFFFF8FFF)==(s2->GetClusterIndex2(i)&0xFFFF8FFF) && s1->GetClusterIndex2(i)>0) {
1989 if ( (s1->GetClusterIndex2(i))==(s2->GetClusterIndex2(i)) && s1->GetClusterIndex2(i)>0) {
1990 AliTPCTrackerPoint *p1 = s1->GetTrackPoint(i);
1991 AliTPCTrackerPoint *p2 = s2->GetTrackPoint(i);;
1992 if (s1->IsActive()&&s2->IsActive()){
1993 p1->fIsShared = kTRUE;
1994 p2->fIsShared = kTRUE;
2001 for (Int_t i=0;i<4;i++){
2002 if (s1->fOverlapLabels[3*i]==0){
2003 s1->fOverlapLabels[3*i] = s2->GetLabel();
2004 s1->fOverlapLabels[3*i+1] = sumshared;
2005 s1->fOverlapLabels[3*i+2] = s2->GetUniqueID();
2009 for (Int_t i=0;i<4;i++){
2010 if (s2->fOverlapLabels[3*i]==0){
2011 s2->fOverlapLabels[3*i] = s1->GetLabel();
2012 s2->fOverlapLabels[3*i+1] = sumshared;
2013 s2->fOverlapLabels[3*i+2] = s1->GetUniqueID();
2021 void AliTPCtrackerMI::SignShared(TObjArray * arr)
2024 //sort trackss according sectors
2026 for (Int_t i=0; i<arr->GetEntriesFast(); i++) {
2027 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2029 //if (pt) RotateToLocal(pt);
2033 arr->Sort(); // sorting according z
2034 arr->Expand(arr->GetEntries());
2037 Int_t nseed=arr->GetEntriesFast();
2038 for (Int_t i=0; i<nseed; i++) {
2039 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2041 for (Int_t j=0;j<=12;j++){
2042 pt->fOverlapLabels[j] =0;
2045 for (Int_t i=0; i<nseed; i++) {
2046 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2048 if (pt->fRemoval>10) continue;
2049 for (Int_t j=i+1; j<nseed; j++){
2050 AliTPCseed *pt2=(AliTPCseed*)arr->UncheckedAt(j);
2052 if (pt2->fRemoval<=10) {
2053 if ( TMath::Abs(pt->fRelativeSector-pt2->fRelativeSector)>0) break;
2060 void AliTPCtrackerMI::RemoveDouble(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t removalindex)
2063 //sort trackss according sectors
2066 Info("RemoveDouble","Number of tracks before double removal- %d\n",arr->GetEntries());
2069 for (Int_t i=0; i<arr->GetEntriesFast(); i++) {
2070 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2075 arr->Sort(); // sorting according z
2076 arr->Expand(arr->GetEntries());
2078 //reset overlap labels
2080 Int_t nseed=arr->GetEntriesFast();
2081 for (Int_t i=0; i<nseed; i++) {
2082 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2085 for (Int_t j=0;j<=12;j++){
2086 pt->fOverlapLabels[j] =0;
2090 //sign shared tracks
2091 for (Int_t i=0; i<nseed; i++) {
2092 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2094 if (pt->fRemoval>10) continue;
2095 Float_t deltac = pt->GetC()*0.1;
2096 for (Int_t j=i+1; j<nseed; j++){
2097 AliTPCseed *pt2=(AliTPCseed*)arr->UncheckedAt(j);
2099 if (pt2->fRemoval<=10) {
2100 if ( TMath::Abs(pt->fRelativeSector-pt2->fRelativeSector)>0) break;
2101 if (TMath::Abs(pt->GetC() -pt2->GetC())>deltac) continue;
2102 if (TMath::Abs(pt->GetTgl()-pt2->GetTgl())>0.05) continue;
2109 // remove highly shared tracks
2110 for (Int_t i=0; i<nseed; i++) {
2111 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2113 if (pt->fRemoval>10) continue;
2116 for (Int_t j=0;j<4;j++){
2117 sumshared = pt->fOverlapLabels[j*3+1];
2119 Float_t factor = factor1;
2120 if (pt->fRemoval>0) factor = factor2;
2121 if (sumshared/pt->GetNumberOfClusters()>factor){
2122 for (Int_t j=0;j<4;j++){
2123 if (pt->fOverlapLabels[3*j]==0) continue;
2124 if (pt->fOverlapLabels[3*j+1]<5) continue;
2125 if (pt->fRemoval==removalindex) continue;
2126 AliTPCseed * pt2 = (AliTPCseed*)arr->UncheckedAt(pt->fOverlapLabels[3*j+2]);
2128 if (pt2->GetSigma2C()<pt->GetSigma2C()){
2129 // pt->fRemoval = removalindex;
2130 delete arr->RemoveAt(i);
2138 Info("RemoveDouble","Number of tracks after double removal- %d\n",arr->GetEntries());
2147 void AliTPCtrackerMI::SortTracks(TObjArray * arr, Int_t mode) const
2150 //sort tracks in array according mode criteria
2151 Int_t nseed = arr->GetEntriesFast();
2152 for (Int_t i=0; i<nseed; i++) {
2153 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2163 void AliTPCtrackerMI::RemoveUsed(TObjArray * arr, Float_t factor1, Float_t factor2, Int_t removalindex)
2166 //Loop over all tracks and remove "overlaps"
2169 Int_t nseed = arr->GetEntriesFast();
2172 for (Int_t i=0; i<nseed; i++) {
2173 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2175 delete arr->RemoveAt(i);
2179 pt->fBSigned = kFALSE;
2183 nseed = arr->GetEntriesFast();
2190 for (Int_t i=0; i<nseed; i++) {
2191 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2195 Int_t found,foundable,shared;
2197 pt->GetClusterStatistic(0,160,found, foundable,shared,kFALSE);
2199 pt->GetClusterStatistic(0,160,found, foundable,shared,kTRUE);
2201 Double_t factor = factor2;
2202 if (pt->fBConstrain) factor = factor1;
2204 if ((Float_t(shared)/Float_t(found))>factor){
2205 pt->Desactivate(removalindex);
2210 for (Int_t i=0; i<160; i++) {
2211 Int_t index=pt->GetClusterIndex2(i);
2212 if (index<0 || index&0x8000 ) continue;
2213 AliTPCclusterMI *c= pt->fClusterPointer[i];
2215 // if (!c->IsUsed(10)) c->Use(10);
2216 //if (pt->IsActive())
2225 Info("RemoveUsed","\n*****\nNumber of good tracks after shared removal\t%d\n",fNtracks);
2229 void AliTPCtrackerMI::UnsignClusters()
2232 // loop over all clusters and unsign them
2235 for (Int_t sec=0;sec<fkNIS;sec++){
2236 for (Int_t row=0;row<fInnerSec->GetNRows();row++){
2237 AliTPCclusterMI *cl = fInnerSec[sec][row].fClusters1;
2238 for (Int_t icl =0;icl< fInnerSec[sec][row].fN1;icl++)
2239 // if (cl[icl].IsUsed(10))
2241 cl = fInnerSec[sec][row].fClusters2;
2242 for (Int_t icl =0;icl< fInnerSec[sec][row].fN2;icl++)
2243 //if (cl[icl].IsUsed(10))
2248 for (Int_t sec=0;sec<fkNOS;sec++){
2249 for (Int_t row=0;row<fOuterSec->GetNRows();row++){
2250 AliTPCclusterMI *cl = fOuterSec[sec][row].fClusters1;
2251 for (Int_t icl =0;icl< fOuterSec[sec][row].fN1;icl++)
2252 //if (cl[icl].IsUsed(10))
2254 cl = fOuterSec[sec][row].fClusters2;
2255 for (Int_t icl =0;icl< fOuterSec[sec][row].fN2;icl++)
2256 //if (cl[icl].IsUsed(10))
2265 void AliTPCtrackerMI::SignClusters(TObjArray * arr, Float_t fnumber, Float_t fdensity)
2268 //sign clusters to be "used"
2270 // snumber and sdensity sign number of sigmas - bellow mean value to be accepted
2271 // loop over "primaries"
2285 Int_t nseed = arr->GetEntriesFast();
2286 for (Int_t i=0; i<nseed; i++) {
2287 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2291 if (!(pt->IsActive())) continue;
2292 Float_t dens = pt->GetNumberOfClusters()/Float_t(pt->fNFoundable);
2293 if ( (dens>0.7) && (pt->GetNumberOfClusters()>70)){
2295 sumdens2+= dens*dens;
2296 sumn += pt->GetNumberOfClusters();
2297 sumn2 += pt->GetNumberOfClusters()*pt->GetNumberOfClusters();
2298 Float_t chi2 = pt->GetChi2()/pt->GetNumberOfClusters();
2301 sumchi2 +=chi2*chi2;
2306 Float_t mdensity = 0.9;
2307 Float_t meann = 130;
2308 Float_t meanchi = 1;
2309 Float_t sdensity = 0.1;
2310 Float_t smeann = 10;
2311 Float_t smeanchi =0.4;
2315 mdensity = sumdens/sum;
2317 meanchi = sumchi/sum;
2319 sdensity = sumdens2/sum-mdensity*mdensity;
2320 sdensity = TMath::Sqrt(sdensity);
2322 smeann = sumn2/sum-meann*meann;
2323 smeann = TMath::Sqrt(smeann);
2325 smeanchi = sumchi2/sum - meanchi*meanchi;
2326 smeanchi = TMath::Sqrt(smeanchi);
2330 //REMOVE SHORT DELTAS or tracks going out of sensitive volume of TPC
2332 for (Int_t i=0; i<nseed; i++) {
2333 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2337 if (pt->fBSigned) continue;
2338 if (pt->fBConstrain) continue;
2339 //if (!(pt->IsActive())) continue;
2341 Int_t found,foundable,shared;
2342 pt->GetClusterStatistic(0,160,found, foundable,shared);
2343 if (shared/float(found)>0.3) {
2344 if (shared/float(found)>0.9 ){
2345 //delete arr->RemoveAt(i);
2350 Bool_t isok =kFALSE;
2351 if ( (pt->fNShared/pt->GetNumberOfClusters()<0.5) &&pt->GetNumberOfClusters()>60)
2353 if ((TMath::Abs(1/pt->GetC())<100.) && (pt->fNShared/pt->GetNumberOfClusters()<0.7))
2355 if (TMath::Abs(pt->GetZ()/pt->GetX())>1.1)
2357 if ( (TMath::Abs(pt->GetSnp()>0.7) && pt->GetD(0,0)>60.))
2361 for (Int_t i=0; i<160; i++) {
2362 Int_t index=pt->GetClusterIndex2(i);
2363 if (index<0) continue;
2364 AliTPCclusterMI *c= pt->fClusterPointer[i];
2366 //if (!(c->IsUsed(10))) c->Use();
2373 Double_t maxchi = meanchi+2.*smeanchi;
2375 for (Int_t i=0; i<nseed; i++) {
2376 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2380 //if (!(pt->IsActive())) continue;
2381 if (pt->fBSigned) continue;
2382 Double_t chi = pt->GetChi2()/pt->GetNumberOfClusters();
2383 if (chi>maxchi) continue;
2386 Float_t dens = pt->GetNumberOfClusters()/Float_t(pt->fNFoundable);
2388 //sign only tracks with enoug big density at the beginning
2390 if ((pt->GetDensityFirst(40)<0.75) && pt->GetNumberOfClusters()<meann) continue;
2393 Double_t mindens = TMath::Max(double(mdensity-sdensity*fdensity*bfactor),0.65);
2394 Double_t minn = TMath::Max(Int_t(meann-fnumber*smeann*bfactor),50);
2396 // if (pt->fBConstrain) mindens = TMath::Max(mdensity-sdensity*fdensity*bfactor,0.65);
2397 if ( (pt->fRemoval==10) && (pt->GetSnp()>0.8)&&(dens>mindens))
2400 if ((dens>mindens && pt->GetNumberOfClusters()>minn) && chi<maxchi ){
2401 //Int_t noc=pt->GetNumberOfClusters();
2402 pt->fBSigned = kTRUE;
2403 for (Int_t i=0; i<160; i++) {
2405 Int_t index=pt->GetClusterIndex2(i);
2406 if (index<0) continue;
2407 AliTPCclusterMI *c= pt->fClusterPointer[i];
2409 // if (!(c->IsUsed(10))) c->Use();
2414 // gLastCheck = nseed;
2422 void AliTPCtrackerMI::StopNotActive(TObjArray * arr, Int_t row0, Float_t th0, Float_t th1, Float_t th2) const
2424 // stop not active tracks
2425 // take th1 as threshold for number of founded to number of foundable on last 10 active rows
2426 // take th2 as threshold for number of founded to number of foundable on last 20 active rows
2427 Int_t nseed = arr->GetEntriesFast();
2429 for (Int_t i=0; i<nseed; i++) {
2430 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
2434 if (!(pt->IsActive())) continue;
2435 StopNotActive(pt,row0,th0, th1,th2);
2441 void AliTPCtrackerMI::StopNotActive(AliTPCseed * seed, Int_t row0, Float_t th0, Float_t th1,
2444 // stop not active tracks
2445 // take th1 as threshold for number of founded to number of foundable on last 10 active rows
2446 // take th2 as threshold for number of founded to number of foundable on last 20 active rows
2449 Int_t foundable = 0;
2450 Int_t maxindex = seed->fLastPoint; //last foundable row
2451 if (seed->fNFoundable*th0 > seed->GetNumberOfClusters()) {
2452 seed->Desactivate(10) ;
2456 for (Int_t i=row0; i<maxindex; i++){
2457 Int_t index = seed->GetClusterIndex2(i);
2458 if (index!=-1) foundable++;
2460 if (foundable<=30) sumgood1++;
2461 if (foundable<=50) {
2468 if (foundable>=30.){
2469 if (sumgood1<(th1*30.)) seed->Desactivate(10);
2472 if (sumgood2<(th2*50.)) seed->Desactivate(10);
2476 Int_t AliTPCtrackerMI::RefitInward(AliESD *event)
2479 // back propagation of ESD tracks
2485 //PrepareForProlongation(fSeeds,1);
2486 PropagateForward2(fSeeds);
2488 Int_t nseed = fSeeds->GetEntriesFast();
2489 for (Int_t i=0;i<nseed;i++){
2490 AliTPCseed * seed = (AliTPCseed*) fSeeds->UncheckedAt(i);
2491 if (!seed) continue;
2492 seed->PropagateTo(fParam->GetInnerRadiusLow());
2493 AliESDtrack *esd=event->GetTrack(i);
2494 seed->CookdEdx(0.02,0.6);
2495 CookLabel(seed,0.1); //For comparison only
2496 if (seed->GetNumberOfClusters()>60){
2497 esd->UpdateTrackParams(seed,AliESDtrack::kTPCrefit);
2501 //printf("problem\n");
2504 Info("RefitInward","Number of refitted tracks %d",ntracks);
2511 Int_t AliTPCtrackerMI::PropagateBack(AliESD *event)
2514 // back propagation of ESD tracks
2520 PropagateBack(fSeeds);
2521 Int_t nseed = fSeeds->GetEntriesFast();
2523 for (Int_t i=0;i<nseed;i++){
2524 AliTPCseed * seed = (AliTPCseed*) fSeeds->UncheckedAt(i);
2525 if (!seed) continue;
2526 AliESDtrack *esd=event->GetTrack(i);
2527 seed->CookdEdx(0.02,0.6);
2528 CookLabel(seed,0.1); //For comparison only
2529 if (seed->GetNumberOfClusters()>60){
2530 esd->UpdateTrackParams(seed,AliESDtrack::kTPCout);
2534 Info("PropagateBack","Number of back propagated tracks %d",ntracks);
2541 void AliTPCtrackerMI::DeleteSeeds()
2545 Int_t nseed = fSeeds->GetEntriesFast();
2546 for (Int_t i=0;i<nseed;i++){
2547 AliTPCseed * seed = (AliTPCseed*)fSeeds->At(i);
2548 if (seed) delete fSeeds->RemoveAt(i);
2554 void AliTPCtrackerMI::ReadSeeds(AliESD *event, Int_t direction)
2557 //read seeds from the event
2559 Int_t nentr=event->GetNumberOfTracks();
2561 Info("ReadSeeds", "Number of ESD tracks: %d\n", nentr);
2566 fSeeds = new TObjArray(nentr);
2570 for (Int_t i=0; i<nentr; i++) {
2571 AliESDtrack *esd=event->GetTrack(i);
2572 ULong_t status=esd->GetStatus();
2573 AliTPCtrack t(*esd);
2574 AliTPCseed *seed = new AliTPCseed(t,t.GetAlpha());
2575 if ((status==AliESDtrack::kTPCin)&&(direction==1)) seed->ResetCovariance();
2576 if ( direction ==2 &&(status & AliESDtrack::kTRDrefit) == 0 ) seed->ResetCovariance();
2577 if ( direction ==2 && ((status & AliESDtrack::kTPCout) == 0) ) {
2582 if ( direction ==2 &&(status & AliESDtrack::kTRDrefit) > 0 ) {
2583 Double_t par0[5],par1[5],x;
2584 esd->GetInnerExternalParameters(x,par0);
2585 esd->GetExternalParameters(x,par1);
2586 Double_t delta1 = TMath::Abs(par0[4]-par1[4])/(0.000000001+TMath::Abs(par0[4]+par1[4]));
2587 Double_t delta2 = TMath::Abs(par0[3]-par1[3]);
2588 //reset covariance if suspicious
2589 if ( (delta1>0.1) || (delta2>0.006))
2590 seed->ResetCovariance();
2595 // rotate to the local coordinate system
2597 fSectors=fInnerSec; fN=fkNIS;
2599 Double_t alpha=seed->GetAlpha() - fSectors->GetAlphaShift();
2600 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
2601 if (alpha < 0. ) alpha += 2.*TMath::Pi();
2602 Int_t ns=Int_t(alpha/fSectors->GetAlpha())%fN;
2603 alpha =ns*fSectors->GetAlpha() + fSectors->GetAlphaShift();
2604 if (alpha<-TMath::Pi()) alpha += 2*TMath::Pi();
2605 if (alpha>=TMath::Pi()) alpha -= 2*TMath::Pi();
2606 alpha-=seed->GetAlpha();
2607 if (!seed->Rotate(alpha)) {
2613 //seed->PropagateTo(fSectors->GetX(0));
2615 // Int_t index = esd->GetTPCindex();
2616 //AliTPCseed * seed2= (AliTPCseed*)fSeeds->At(index);
2617 //if (direction==2){
2618 // AliTPCseed * seed2 = ReSeed(seed,0.,0.5,1.);
2626 for (Int_t irow=0;irow<160;irow++){
2627 Int_t index = seed->GetClusterIndex2(irow);
2630 AliTPCclusterMI * cl = GetClusterMI(index);
2631 seed->fClusterPointer[irow] = cl;
2633 if ((index & 0x8000)==0){
2634 cl->Use(10); // accepted cluster
2636 cl->Use(6); // close cluster not accepted
2639 Info("ReadSeeds","Not found cluster");
2643 fSeeds->AddAt(seed,i);
2649 //_____________________________________________________________________________
2650 void AliTPCtrackerMI::MakeSeeds3(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4],
2651 Float_t deltay, Int_t ddsec) {
2652 //-----------------------------------------------------------------
2653 // This function creates track seeds.
2654 // SEEDING WITH VERTEX CONSTRAIN
2655 //-----------------------------------------------------------------
2656 // cuts[0] - fP4 cut
2657 // cuts[1] - tan(phi) cut
2658 // cuts[2] - zvertex cut
2659 // cuts[3] - fP3 cut
2667 Double_t x[5], c[15];
2668 // Int_t di = i1-i2;
2670 AliTPCseed * seed = new AliTPCseed;
2671 Double_t alpha=fSectors->GetAlpha(), shift=fSectors->GetAlphaShift();
2672 Double_t cs=cos(alpha), sn=sin(alpha);
2674 // Double_t x1 =fOuterSec->GetX(i1);
2675 //Double_t xx2=fOuterSec->GetX(i2);
2677 Double_t x1 =GetXrow(i1);
2678 Double_t xx2=GetXrow(i2);
2680 Double_t x3=GetX(), y3=GetY(), z3=GetZ();
2682 Int_t imiddle = (i2+i1)/2; //middle pad row index
2683 Double_t xm = GetXrow(imiddle); // radius of middle pad-row
2684 const AliTPCRow& krm=GetRow(sec,imiddle); //middle pad -row
2688 const AliTPCRow& kr1=GetRow(ns,i1);
2689 Double_t ymax = GetMaxY(i1)-kr1.fDeadZone-1.5;
2690 Double_t ymaxm = GetMaxY(imiddle)-kr1.fDeadZone-1.5;
2693 // change cut on curvature if it can't reach this layer
2694 // maximal curvature set to reach it
2695 Double_t dvertexmax = TMath::Sqrt((x1-x3)*(x1-x3)+(ymax+5-y3)*(ymax+5-y3));
2696 if (dvertexmax*0.5*cuts[0]>0.85){
2697 cuts[0] = 0.85/(dvertexmax*0.5+1.);
2699 Double_t r2min = 1/(cuts[0]*cuts[0]); //minimal square of radius given by cut
2702 if (deltay>0) ddsec = 0;
2703 // loop over clusters
2704 for (Int_t is=0; is < kr1; is++) {
2706 if (kr1[is]->IsUsed(10)) continue;
2707 Double_t y1=kr1[is]->GetY(), z1=kr1[is]->GetZ();
2708 //if (TMath::Abs(y1)>ymax) continue;
2710 if (deltay>0 && TMath::Abs(ymax-TMath::Abs(y1))> deltay ) continue; // seed only at the edge
2712 // find possible directions
2713 Float_t anglez = (z1-z3)/(x1-x3);
2714 Float_t extraz = z1 - anglez*(x1-xx2); // extrapolated z
2717 //find rotation angles relative to line given by vertex and point 1
2718 Double_t dvertex2 = (x1-x3)*(x1-x3)+(y1-y3)*(y1-y3);
2719 Double_t dvertex = TMath::Sqrt(dvertex2);
2720 Double_t angle13 = TMath::ATan((y1-y3)/(x1-x3));
2721 Double_t cs13 = cos(-angle13), sn13 = sin(-angle13);
2724 // loop over 2 sectors
2730 Double_t dddz1=0; // direction of delta inclination in z axis
2737 for (Int_t dsec = dsec1; dsec<=dsec2;dsec++){
2738 Int_t sec2 = sec + dsec;
2740 // AliTPCRow& kr2 = fOuterSec[(sec2+fkNOS)%fkNOS][i2];
2741 //AliTPCRow& kr2m = fOuterSec[(sec2+fkNOS)%fkNOS][imiddle];
2742 AliTPCRow& kr2 = GetRow((sec2+fkNOS)%fkNOS,i2);
2743 AliTPCRow& kr2m = GetRow((sec2+fkNOS)%fkNOS,imiddle);
2744 Int_t index1 = TMath::Max(kr2.Find(extraz-0.6-dddz1*TMath::Abs(z1)*0.05)-1,0);
2745 Int_t index2 = TMath::Min(kr2.Find(extraz+0.6+dddz2*TMath::Abs(z1)*0.05)+1,kr2);
2747 // rotation angles to p1-p3
2748 Double_t cs13r = cos(-angle13+dsec*alpha)/dvertex, sn13r = sin(-angle13+dsec*alpha)/dvertex;
2749 Double_t x2, y2, z2;
2751 // Double_t dymax = maxangle*TMath::Abs(x1-xx2);
2754 Double_t dxx0 = (xx2-x3)*cs13r;
2755 Double_t dyy0 = (xx2-x3)*sn13r;
2756 for (Int_t js=index1; js < index2; js++) {
2757 const AliTPCclusterMI *kcl = kr2[js];
2758 if (kcl->IsUsed(10)) continue;
2760 //calcutate parameters
2762 Double_t yy0 = dyy0 +(kcl->GetY()-y3)*cs13r;
2764 if (TMath::Abs(yy0)<0.000001) continue;
2765 Double_t xx0 = dxx0 -(kcl->GetY()-y3)*sn13r;
2766 Double_t y0 = 0.5*(xx0*xx0+yy0*yy0-xx0)/yy0;
2767 Double_t r02 = (0.25+y0*y0)*dvertex2;
2768 //curvature (radius) cut
2769 if (r02<r2min) continue;
2773 Double_t c0 = 1/TMath::Sqrt(r02);
2777 //Double_t dfi0 = 2.*TMath::ASin(dvertex*c0*0.5);
2778 //Double_t dfi1 = 2.*TMath::ASin(TMath::Sqrt(yy0*yy0+(1-xx0)*(1-xx0))*dvertex*c0*0.5);
2779 Double_t dfi0 = 2.*AliTPCFastMath::FastAsin(dvertex*c0*0.5);
2780 Double_t dfi1 = 2.*AliTPCFastMath::FastAsin(TMath::Sqrt(yy0*yy0+(1-xx0)*(1-xx0))*dvertex*c0*0.5);
2783 Double_t z0 = kcl->GetZ();
2784 Double_t zzzz2 = z1-(z1-z3)*dfi1/dfi0;
2785 if (TMath::Abs(zzzz2-z0)>0.5) continue;
2788 Double_t dip = (z1-z0)*c0/dfi1;
2789 Double_t x0 = (0.5*cs13+y0*sn13)*dvertex*c0;
2800 x2= xx2*cs-y2*sn*dsec;
2801 y2=+xx2*sn*dsec+y2*cs;
2811 // do we have cluster at the middle ?
2813 GetProlongation(x1,xm,x,ym,zm);
2815 AliTPCclusterMI * cm=0;
2816 if (TMath::Abs(ym)-ymaxm<0){
2817 cm = krm.FindNearest2(ym,zm,1.0,0.6,dummy);
2818 if ((!cm) || (cm->IsUsed(10))) {
2823 // rotate y1 to system 0
2824 // get state vector in rotated system
2825 Double_t yr1 = (-0.5*sn13+y0*cs13)*dvertex*c0;
2826 Double_t xr2 = x0*cs+yr1*sn*dsec;
2827 Double_t xr[5]={kcl->GetY(),kcl->GetZ(), xr2, dip, c0};
2829 GetProlongation(xx2,xm,xr,ym,zm);
2830 if (TMath::Abs(ym)-ymaxm<0){
2831 cm = kr2m.FindNearest2(ym,zm,1.0,0.6,dummy);
2832 if ((!cm) || (cm->IsUsed(10))) {
2842 dym = ym - cm->GetY();
2843 dzm = zm - cm->GetZ();
2850 Double_t sy1=kr1[is]->GetSigmaY2()*2., sz1=kr1[is]->GetSigmaZ2()*2.;
2851 Double_t sy2=kcl->GetSigmaY2()*2., sz2=kcl->GetSigmaZ2()*2.;
2852 //Double_t sy3=400*3./12., sy=0.1, sz=0.1;
2853 Double_t sy3=25000*x[4]*x[4]+0.1, sy=0.1, sz=0.1;
2854 //Double_t sy3=25000*x[4]*x[4]*60+0.5, sy=0.1, sz=0.1;
2856 Double_t f40=(F1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
2857 Double_t f42=(F1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
2858 Double_t f43=(F1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
2859 Double_t f20=(F2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
2860 Double_t f22=(F2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
2861 Double_t f23=(F2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
2863 Double_t f30=(F3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
2864 Double_t f31=(F3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
2865 Double_t f32=(F3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
2866 Double_t f34=(F3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
2870 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
2871 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
2872 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
2873 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
2874 c[13]=f30*sy1*f40+f32*sy2*f42;
2875 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
2877 // if (!BuildSeed(kr1[is],kcl,0,x1,x2,x3,x,c)) continue;
2879 UInt_t index=kr1.GetIndex(is);
2880 AliTPCseed *track=new(seed) AliTPCseed(index, x, c, x1, ns*alpha+shift);
2882 track->fIsSeeding = kTRUE;
2889 FollowProlongation(*track, (i1+i2)/2,1);
2890 Int_t foundable,found,shared;
2891 track->GetClusterStatistic((i1+i2)/2,i1, found, foundable, shared, kTRUE);
2892 if ((found<0.55*foundable) || shared>0.5*found || (track->GetSigmaY2()+track->GetSigmaZ2())>0.5){
2894 seed->~AliTPCseed();
2900 FollowProlongation(*track, i2,1);
2904 track->fBConstrain =1;
2905 // track->fLastPoint = i1+fInnerSec->GetNRows(); // first cluster in track position
2906 track->fLastPoint = i1; // first cluster in track position
2907 track->fFirstPoint = track->fLastPoint;
2909 if (track->GetNumberOfClusters()<(i1-i2)*0.5 ||
2910 track->GetNumberOfClusters() < track->fNFoundable*0.6 ||
2911 track->fNShared>0.4*track->GetNumberOfClusters() ) {
2913 seed->~AliTPCseed();
2917 // Z VERTEX CONDITION
2919 zv = track->GetZ()+track->GetTgl()/track->GetC()*
2920 ( asin(-track->GetEta()) - asin(track->GetX()*track->GetC()-track->GetEta()));
2921 if (TMath::Abs(zv-z3)>cuts[2]) {
2922 FollowProlongation(*track, TMath::Max(i2-20,0));
2923 zv = track->GetZ()+track->GetTgl()/track->GetC()*
2924 ( asin(-track->GetEta()) - asin(track->GetX()*track->GetC()-track->GetEta()));
2925 if (TMath::Abs(zv-z3)>cuts[2]){
2926 FollowProlongation(*track, TMath::Max(i2-40,0));
2927 zv = track->GetZ()+track->GetTgl()/track->GetC()*
2928 ( asin(-track->GetEta()) - asin(track->GetX()*track->GetC()-track->GetEta()));
2929 if (TMath::Abs(zv-z3)>cuts[2] &&(track->GetNumberOfClusters() > track->fNFoundable*0.7)){
2930 // make seed without constrain
2931 AliTPCseed * track2 = MakeSeed(track,0.2,0.5,1.);
2932 FollowProlongation(*track2, i2,1);
2933 track2->fBConstrain = kFALSE;
2934 track2->fSeedType = 1;
2935 arr->AddLast(track2);
2937 seed->~AliTPCseed();
2942 seed->~AliTPCseed();
2949 track->fSeedType =0;
2950 arr->AddLast(track);
2951 seed = new AliTPCseed;
2953 // don't consider other combinations
2954 if (track->GetNumberOfClusters() > track->fNFoundable*0.8)
2960 Info("MakeSeeds3","\nSeeding statistic:\t%d\t%d\t%d\t%d\t%d\t%d",nin0,nin1,nin2,nin,nout1,nout2);
2966 void AliTPCtrackerMI::MakeSeeds5(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t cuts[4],
2971 //-----------------------------------------------------------------
2972 // This function creates track seeds.
2973 //-----------------------------------------------------------------
2974 // cuts[0] - fP4 cut
2975 // cuts[1] - tan(phi) cut
2976 // cuts[2] - zvertex cut
2977 // cuts[3] - fP3 cut
2987 Double_t x[5], c[15];
2989 // make temporary seed
2990 AliTPCseed * seed = new AliTPCseed;
2991 Double_t alpha=fOuterSec->GetAlpha(), shift=fOuterSec->GetAlphaShift();
2992 // Double_t cs=cos(alpha), sn=sin(alpha);
2997 Double_t x1 = GetXrow(i1-1);
2998 const AliTPCRow& kr1=GetRow(sec,i1-1);
2999 Double_t y1max = GetMaxY(i1-1)-kr1.fDeadZone-1.5;
3001 Double_t x1p = GetXrow(i1);
3002 const AliTPCRow& kr1p=GetRow(sec,i1);
3004 Double_t x1m = GetXrow(i1-2);
3005 const AliTPCRow& kr1m=GetRow(sec,i1-2);
3008 //last 3 padrow for seeding
3009 AliTPCRow& kr3 = GetRow((sec+fkNOS)%fkNOS,i1-7);
3010 Double_t x3 = GetXrow(i1-7);
3011 // Double_t y3max= GetMaxY(i1-7)-kr3.fDeadZone-1.5;
3013 AliTPCRow& kr3p = GetRow((sec+fkNOS)%fkNOS,i1-6);
3014 Double_t x3p = GetXrow(i1-6);
3016 AliTPCRow& kr3m = GetRow((sec+fkNOS)%fkNOS,i1-8);
3017 Double_t x3m = GetXrow(i1-8);
3022 Int_t im = i1-4; //middle pad row index
3023 Double_t xm = GetXrow(im); // radius of middle pad-row
3024 const AliTPCRow& krm=GetRow(sec,im); //middle pad -row
3025 // Double_t ymmax = GetMaxY(im)-kr1.fDeadZone-1.5;
3028 Double_t deltax = x1-x3;
3029 Double_t dymax = deltax*cuts[1];
3030 Double_t dzmax = deltax*cuts[3];
3032 // loop over clusters
3033 for (Int_t is=0; is < kr1; is++) {
3035 if (kr1[is]->IsUsed(10)) continue;
3036 Double_t y1=kr1[is]->GetY(), z1=kr1[is]->GetZ();
3038 if (deltay>0 && TMath::Abs(y1max-TMath::Abs(y1))> deltay ) continue; // seed only at the edge
3040 Int_t index1 = TMath::Max(kr3.Find(z1-dzmax)-1,0);
3041 Int_t index2 = TMath::Min(kr3.Find(z1+dzmax)+1,kr3);
3047 for (Int_t js=index1; js < index2; js++) {
3048 const AliTPCclusterMI *kcl = kr3[js];
3049 if (kcl->IsUsed(10)) continue;
3051 // apply angular cuts
3052 if (TMath::Abs(y1-y3)>dymax) continue;
3055 if (TMath::Abs(z1-z3)>dzmax) continue;
3057 Double_t angley = (y1-y3)/(x1-x3);
3058 Double_t anglez = (z1-z3)/(x1-x3);
3060 Double_t erry = TMath::Abs(angley)*(x1-x1m)*0.5+0.5;
3061 Double_t errz = TMath::Abs(anglez)*(x1-x1m)*0.5+0.5;
3063 Double_t yyym = angley*(xm-x1)+y1;
3064 Double_t zzzm = anglez*(xm-x1)+z1;
3066 const AliTPCclusterMI *kcm = krm.FindNearest2(yyym,zzzm,erry,errz,index);
3068 if (kcm->IsUsed(10)) continue;
3070 erry = TMath::Abs(angley)*(x1-x1m)*0.4+0.5;
3071 errz = TMath::Abs(anglez)*(x1-x1m)*0.4+0.5;
3078 // look around first
3079 const AliTPCclusterMI *kc1m = kr1m.FindNearest2(angley*(x1m-x1)+y1,
3085 if (kc1m->IsUsed(10)) used++;
3087 const AliTPCclusterMI *kc1p = kr1p.FindNearest2(angley*(x1p-x1)+y1,
3093 if (kc1p->IsUsed(10)) used++;
3095 if (used>1) continue;
3096 if (found<1) continue;
3100 const AliTPCclusterMI *kc3m = kr3m.FindNearest2(angley*(x3m-x3)+y3,
3106 if (kc3m->IsUsed(10)) used++;
3110 const AliTPCclusterMI *kc3p = kr3p.FindNearest2(angley*(x3p-x3)+y3,
3116 if (kc3p->IsUsed(10)) used++;
3120 if (used>1) continue;
3121 if (found<3) continue;
3131 x[4]=F1(x1,y1,x2,y2,x3,y3);
3132 //if (TMath::Abs(x[4]) >= cuts[0]) continue;
3135 x[2]=F2(x1,y1,x2,y2,x3,y3);
3138 x[3]=F3n(x1,y1,x2,y2,z1,z2,x[4]);
3139 //if (TMath::Abs(x[3]) > cuts[3]) continue;
3143 Double_t sy1=0.1, sz1=0.1;
3144 Double_t sy2=0.1, sz2=0.1;
3145 Double_t sy3=0.1, sy=0.1, sz=0.1;
3147 Double_t f40=(F1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
3148 Double_t f42=(F1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
3149 Double_t f43=(F1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
3150 Double_t f20=(F2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
3151 Double_t f22=(F2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
3152 Double_t f23=(F2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
3154 Double_t f30=(F3(x1,y1+sy,x2,y2,z1,z2)-x[3])/sy;
3155 Double_t f31=(F3(x1,y1,x2,y2,z1+sz,z2)-x[3])/sz;
3156 Double_t f32=(F3(x1,y1,x2,y2+sy,z1,z2)-x[3])/sy;
3157 Double_t f34=(F3(x1,y1,x2,y2,z1,z2+sz)-x[3])/sz;
3161 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
3162 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
3163 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
3164 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
3165 c[13]=f30*sy1*f40+f32*sy2*f42;
3166 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
3168 // if (!BuildSeed(kr1[is],kcl,0,x1,x2,x3,x,c)) continue;
3170 UInt_t index=kr1.GetIndex(is);
3171 AliTPCseed *track=new(seed) AliTPCseed(index, x, c, x1, sec*alpha+shift);
3173 track->fIsSeeding = kTRUE;
3176 FollowProlongation(*track, i1-7,1);
3177 if (track->GetNumberOfClusters() < track->fNFoundable*0.75 ||
3178 track->fNShared>0.6*track->GetNumberOfClusters() || ( track->GetSigmaY2()+ track->GetSigmaZ2())>0.6){
3180 seed->~AliTPCseed();
3186 FollowProlongation(*track, i2,1);
3187 track->fBConstrain =0;
3188 track->fLastPoint = i1+fInnerSec->GetNRows(); // first cluster in track position
3189 track->fFirstPoint = track->fLastPoint;
3191 if (track->GetNumberOfClusters()<(i1-i2)*0.5 ||
3192 track->GetNumberOfClusters()<track->fNFoundable*0.7 ||
3193 track->fNShared>2. || track->GetChi2()/track->GetNumberOfClusters()>6 || ( track->GetSigmaY2()+ track->GetSigmaZ2())>0.5 ) {
3195 seed->~AliTPCseed();
3200 FollowProlongation(*track, TMath::Max(i2-10,0),1);
3201 AliTPCseed * track2 = MakeSeed(track,0.2,0.5,0.9);
3202 FollowProlongation(*track2, i2,1);
3203 track2->fBConstrain = kFALSE;
3204 track2->fSeedType = 4;
3205 arr->AddLast(track2);
3207 seed->~AliTPCseed();
3211 //arr->AddLast(track);
3212 //seed = new AliTPCseed;
3218 Info("MakeSeeds5","\nSeeding statiistic:\t%d\t%d\t%d\t%d\t%d\t%d",nin0,nin1,nin2,nin,nout1,nout2,nout3);
3224 //_____________________________________________________________________________
3225 void AliTPCtrackerMI::MakeSeeds2(TObjArray * arr, Int_t sec, Int_t i1, Int_t i2, Float_t */*cuts[4]*/,
3226 Float_t deltay, Bool_t /*bconstrain*/) {
3227 //-----------------------------------------------------------------
3228 // This function creates track seeds - without vertex constraint
3229 //-----------------------------------------------------------------
3230 // cuts[0] - fP4 cut - not applied
3231 // cuts[1] - tan(phi) cut
3232 // cuts[2] - zvertex cut - not applied
3233 // cuts[3] - fP3 cut
3243 Double_t alpha=fOuterSec->GetAlpha(), shift=fOuterSec->GetAlphaShift();
3244 // Double_t cs=cos(alpha), sn=sin(alpha);
3245 Int_t row0 = (i1+i2)/2;
3246 Int_t drow = (i1-i2)/2;
3247 const AliTPCRow& kr0=fSectors[sec][row0];
3250 AliTPCpolyTrack polytrack;
3251 Int_t nclusters=fSectors[sec][row0];
3252 AliTPCseed * seed = new AliTPCseed;
3257 for (Int_t is=0; is < nclusters; is++) { //LOOP over clusters
3259 Int_t nfoundable =0;
3260 for (Int_t iter =1; iter<2; iter++){ //iterations
3261 const AliTPCRow& krm=fSectors[sec][row0-iter];
3262 const AliTPCRow& krp=fSectors[sec][row0+iter];
3263 const AliTPCclusterMI * cl= kr0[is];
3265 if (cl->IsUsed(10)) {
3271 Double_t x = kr0.GetX();
3272 // Initialization of the polytrack
3277 Double_t y0= cl->GetY();
3278 Double_t z0= cl->GetZ();
3282 Double_t ymax = fSectors->GetMaxY(row0)-kr0.fDeadZone-1.5;
3283 if (deltay>0 && TMath::Abs(ymax-TMath::Abs(y0))> deltay ) continue; // seed only at the edge
3285 erry = (0.5)*cl->GetSigmaY2()/TMath::Sqrt(cl->GetQ())*6;
3286 errz = (0.5)*cl->GetSigmaZ2()/TMath::Sqrt(cl->GetQ())*6;
3287 polytrack.AddPoint(x,y0,z0,erry, errz);
3290 if (cl->IsUsed(10)) sumused++;
3293 Float_t roady = (5*TMath::Sqrt(cl->GetSigmaY2()+0.2)+1.)*iter;
3294 Float_t roadz = (5*TMath::Sqrt(cl->GetSigmaZ2()+0.2)+1.)*iter;
3297 AliTPCclusterMI * cl1 = krm.FindNearest(y0,z0,roady,roadz);
3298 if (cl1 && TMath::Abs(ymax-TMath::Abs(y0))) {
3299 erry = (0.5)*cl1->GetSigmaY2()/TMath::Sqrt(cl1->GetQ())*3;
3300 errz = (0.5)*cl1->GetSigmaZ2()/TMath::Sqrt(cl1->GetQ())*3;
3301 if (cl1->IsUsed(10)) sumused++;
3302 polytrack.AddPoint(x,cl1->GetY(),cl1->GetZ(),erry,errz);
3306 AliTPCclusterMI * cl2 = krp.FindNearest(y0,z0,roady,roadz);
3308 erry = (0.5)*cl2->GetSigmaY2()/TMath::Sqrt(cl2->GetQ())*3;
3309 errz = (0.5)*cl2->GetSigmaZ2()/TMath::Sqrt(cl2->GetQ())*3;
3310 if (cl2->IsUsed(10)) sumused++;
3311 polytrack.AddPoint(x,cl2->GetY(),cl2->GetZ(),erry,errz);
3314 if (sumused>0) continue;
3316 polytrack.UpdateParameters();
3322 nfoundable = polytrack.GetN();
3323 nfound = nfoundable;
3325 for (Int_t ddrow = iter+1; ddrow<drow;ddrow++){
3326 Float_t maxdist = 0.8*(1.+3./(ddrow));
3327 for (Int_t delta = -1;delta<=1;delta+=2){
3328 Int_t row = row0+ddrow*delta;
3329 kr = &(fSectors[sec][row]);
3330 Double_t xn = kr->GetX();
3331 Double_t ymax = fSectors->GetMaxY(row)-kr->fDeadZone-1.5;
3332 polytrack.GetFitPoint(xn,yn,zn);
3333 if (TMath::Abs(yn)>ymax) continue;
3335 AliTPCclusterMI * cln = kr->FindNearest(yn,zn,roady,roadz);
3337 Float_t dist = TMath::Sqrt( (yn-cln->GetY())*(yn-cln->GetY())+(zn-cln->GetZ())*(zn-cln->GetZ()));
3340 erry = (dist+0.3)*cln->GetSigmaY2()/TMath::Sqrt(cln->GetQ())*(1.+1./(ddrow));
3341 errz = (dist+0.3)*cln->GetSigmaZ2()/TMath::Sqrt(cln->GetQ())*(1.+1./(ddrow));
3342 if (cln->IsUsed(10)) {
3343 // printf("used\n");
3351 polytrack.AddPoint(xn,cln->GetY(),cln->GetZ(),erry, errz);
3356 if ( (sumused>3) || (sumused>0.5*nfound) || (nfound<0.6*nfoundable)) break;
3357 polytrack.UpdateParameters();
3360 if ( (sumused>3) || (sumused>0.5*nfound)) {
3361 //printf("sumused %d\n",sumused);
3366 polytrack.GetFitDerivation(kr0.GetX(),dy,dz);
3367 AliTPCpolyTrack track2;
3369 polytrack.Refit(track2,0.5+TMath::Abs(dy)*0.3,0.4+TMath::Abs(dz)*0.3);
3370 if (track2.GetN()<0.5*nfoundable) continue;
3373 if ((nfound>0.6*nfoundable) &&( nfoundable>0.4*(i1-i2))) {
3375 // test seed with and without constrain
3376 for (Int_t constrain=0; constrain<=0;constrain++){
3377 // add polytrack candidate
3379 Double_t x[5], c[15];
3380 Double_t x1,x2,x3,y1,y2,y3,z1,z2,z3;
3381 track2.GetBoundaries(x3,x1);
3383 track2.GetFitPoint(x1,y1,z1);
3384 track2.GetFitPoint(x2,y2,z2);
3385 track2.GetFitPoint(x3,y3,z3);
3387 //is track pointing to the vertex ?
3390 polytrack.GetFitPoint(x0,y0,z0);
3403 x[4]=F1(x1,y1,x2,y2,x3,y3);
3405 // if (TMath::Abs(x[4]) >= cuts[0]) continue; //
3406 x[2]=F2(x1,y1,x2,y2,x3,y3);
3408 //if (TMath::Abs(x[4]*x1-x[2]) >= cuts[1]) continue;
3409 //x[3]=F3(x1,y1,x2,y2,z1,z2);
3410 x[3]=F3n(x1,y1,x3,y3,z1,z3,x[4]);
3411 //if (TMath::Abs(x[3]) > cuts[3]) continue;
3414 Double_t sy =0.1, sz =0.1;
3415 Double_t sy1=0.02, sz1=0.02;
3416 Double_t sy2=0.02, sz2=0.02;
3420 sy3=25000*x[4]*x[4]+0.1, sy=0.1, sz=0.1;
3423 Double_t f40=(F1(x1,y1+sy,x2,y2,x3,y3)-x[4])/sy;
3424 Double_t f42=(F1(x1,y1,x2,y2+sy,x3,y3)-x[4])/sy;
3425 Double_t f43=(F1(x1,y1,x2,y2,x3,y3+sy)-x[4])/sy;
3426 Double_t f20=(F2(x1,y1+sy,x2,y2,x3,y3)-x[2])/sy;
3427 Double_t f22=(F2(x1,y1,x2,y2+sy,x3,y3)-x[2])/sy;
3428 Double_t f23=(F2(x1,y1,x2,y2,x3,y3+sy)-x[2])/sy;
3430 Double_t f30=(F3(x1,y1+sy,x3,y3,z1,z3)-x[3])/sy;
3431 Double_t f31=(F3(x1,y1,x3,y3,z1+sz,z3)-x[3])/sz;
3432 Double_t f32=(F3(x1,y1,x3,y3+sy,z1,z3)-x[3])/sy;
3433 Double_t f34=(F3(x1,y1,x3,y3,z1,z3+sz)-x[3])/sz;
3438 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
3439 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
3440 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
3441 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
3442 c[13]=f30*sy1*f40+f32*sy2*f42;
3443 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
3445 //Int_t row1 = fSectors->GetRowNumber(x1);
3446 Int_t row1 = GetRowNumber(x1);
3450 AliTPCseed *track=new (seed) AliTPCseed(index, x, c, x1, sec*alpha+shift);
3451 track->fIsSeeding = kTRUE;
3452 Int_t rc=FollowProlongation(*track, i2);
3453 if (constrain) track->fBConstrain =1;
3455 track->fBConstrain =0;
3456 track->fLastPoint = row1+fInnerSec->GetNRows(); // first cluster in track position
3457 track->fFirstPoint = track->fLastPoint;
3459 if (rc==0 || track->GetNumberOfClusters()<(i1-i2)*0.5 ||
3460 track->GetNumberOfClusters() < track->fNFoundable*0.6 ||
3461 track->fNShared>0.4*track->GetNumberOfClusters()) {
3464 seed->~AliTPCseed();
3467 arr->AddLast(track);
3468 seed = new AliTPCseed;
3472 } // if accepted seed
3475 Info("MakeSeeds2","\nSeeding statiistic:\t%d\t%d\t%d\t%d",nin0,nin1,nin2,nin3);
3481 AliTPCseed *AliTPCtrackerMI::MakeSeed(AliTPCseed *track, Float_t r0, Float_t r1, Float_t r2)
3485 //reseed using track points
3486 Int_t p0 = int(r0*track->GetNumberOfClusters()); // point 0
3487 Int_t p1 = int(r1*track->GetNumberOfClusters());
3488 Int_t p2 = int(r2*track->GetNumberOfClusters()); // last point
3490 Double_t x0[3],x1[3],x2[3];
3495 // find track position at given ratio of the length
3496 Int_t sec0, sec1, sec2;
3502 for (Int_t i=0;i<160;i++){
3503 if (track->fClusterPointer[i]){
3505 AliTPCTrackerPoint *trpoint =track->GetTrackPoint(i);
3506 if ( (index<p0) || x0[0]<0 ){
3507 if (trpoint->GetX()>1){
3508 clindex = track->GetClusterIndex2(i);
3510 x0[0] = trpoint->GetX();
3511 x0[1] = trpoint->GetY();
3512 x0[2] = trpoint->GetZ();
3513 sec0 = ((clindex&0xff000000)>>24)%18;
3518 if ( (index<p1) &&(trpoint->GetX()>1)){
3519 clindex = track->GetClusterIndex2(i);
3521 x1[0] = trpoint->GetX();
3522 x1[1] = trpoint->GetY();
3523 x1[2] = trpoint->GetZ();
3524 sec1 = ((clindex&0xff000000)>>24)%18;
3527 if ( (index<p2) &&(trpoint->GetX()>1)){
3528 clindex = track->GetClusterIndex2(i);
3530 x2[0] = trpoint->GetX();
3531 x2[1] = trpoint->GetY();
3532 x2[2] = trpoint->GetZ();
3533 sec2 = ((clindex&0xff000000)>>24)%18;
3540 Double_t alpha, cs,sn, xx2,yy2;
3542 alpha = (sec1-sec2)*fSectors->GetAlpha();
3543 cs = TMath::Cos(alpha);
3544 sn = TMath::Sin(alpha);
3545 xx2= x1[0]*cs-x1[1]*sn;
3546 yy2= x1[0]*sn+x1[1]*cs;
3550 alpha = (sec0-sec2)*fSectors->GetAlpha();
3551 cs = TMath::Cos(alpha);
3552 sn = TMath::Sin(alpha);
3553 xx2= x0[0]*cs-x0[1]*sn;
3554 yy2= x0[0]*sn+x0[1]*cs;
3560 Double_t x[5],c[15];
3564 x[4]=F1(x2[0],x2[1],x1[0],x1[1],x0[0],x0[1]);
3565 // if (x[4]>1) return 0;
3566 x[2]=F2(x2[0],x2[1],x1[0],x1[1],x0[0],x0[1]);
3567 x[3]=F3n(x2[0],x2[1],x0[0],x0[1],x2[2],x0[2],x[4]);
3568 //if (TMath::Abs(x[3]) > 2.2) return 0;
3569 //if (TMath::Abs(x[2]) > 1.99) return 0;
3571 Double_t sy =0.1, sz =0.1;
3573 Double_t sy1=0.02+track->GetSigmaY2(), sz1=0.02+track->GetSigmaZ2();
3574 Double_t sy2=0.01+track->GetSigmaY2(), sz2=0.01+track->GetSigmaZ2();
3575 Double_t sy3=0.01+track->GetSigmaY2();
3577 Double_t f40=(F1(x2[0],x2[1]+sy,x1[0],x1[1],x0[0],x0[1])-x[4])/sy;
3578 Double_t f42=(F1(x2[0],x2[1],x1[0],x1[1]+sy,x0[0],x0[1])-x[4])/sy;
3579 Double_t f43=(F1(x2[0],x2[1],x1[0],x1[1],x0[0],x0[1]+sy)-x[4])/sy;
3580 Double_t f20=(F2(x2[0],x2[1]+sy,x1[0],x1[1],x0[0],x0[1])-x[2])/sy;
3581 Double_t f22=(F2(x2[0],x2[1],x1[0],x1[1]+sy,x0[0],x0[1])-x[2])/sy;
3582 Double_t f23=(F2(x2[0],x2[1],x1[0],x1[1],x0[0],x0[1]+sy)-x[2])/sy;
3584 Double_t f30=(F3(x2[0],x2[1]+sy,x0[0],x0[1],x2[2],x0[2])-x[3])/sy;
3585 Double_t f31=(F3(x2[0],x2[1],x0[0],x0[1],x2[2]+sz,x0[2])-x[3])/sz;
3586 Double_t f32=(F3(x2[0],x2[1],x0[0],x0[1]+sy,x2[2],x0[2])-x[3])/sy;
3587 Double_t f34=(F3(x2[0],x2[1],x0[0],x0[1],x2[2],x0[2]+sz)-x[3])/sz;
3592 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
3593 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
3594 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
3595 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
3596 c[13]=f30*sy1*f40+f32*sy2*f42;
3597 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
3599 // Int_t row1 = fSectors->GetRowNumber(x2[0]);
3600 AliTPCseed *seed=new AliTPCseed(0, x, c, x2[0], sec2*fSectors->GetAlpha()+fSectors->GetAlphaShift());
3601 // Double_t y0,z0,y1,z1, y2,z2;
3602 //seed->GetProlongation(x0[0],y0,z0);
3603 // seed->GetProlongation(x1[0],y1,z1);
3604 //seed->GetProlongation(x2[0],y2,z2);
3606 seed->fLastPoint = pp2;
3607 seed->fFirstPoint = pp2;
3614 AliTPCseed *AliTPCtrackerMI::ReSeed(AliTPCseed *track, Float_t r0, Float_t r1, Float_t r2)
3618 //reseed using founded clusters
3620 // Find the number of clusters
3621 Int_t nclusters = 0;
3622 for (Int_t irow=0;irow<160;irow++){
3623 if (track->GetClusterIndex(irow)>0) nclusters++;
3627 ipos[0] = TMath::Max(int(r0*nclusters),0); // point 0 cluster
3628 ipos[1] = TMath::Min(int(r1*nclusters),nclusters-1); //
3629 ipos[2] = TMath::Min(int(r2*nclusters),nclusters-1); // last point
3633 Int_t row[3],sec[3]={0,0,0};
3635 // find track row position at given ratio of the length
3637 for (Int_t irow=0;irow<160;irow++){
3638 if (track->GetClusterIndex2(irow)<0) continue;
3640 for (Int_t ipoint=0;ipoint<3;ipoint++){
3641 if (index<=ipos[ipoint]) row[ipoint] = irow;
3645 //Get cluster and sector position
3646 for (Int_t ipoint=0;ipoint<3;ipoint++){
3647 Int_t clindex = track->GetClusterIndex2(row[ipoint]);
3648 AliTPCclusterMI * cl = GetClusterMI(clindex);
3651 // AliTPCclusterMI * cl = GetClusterMI(clindex);
3654 sec[ipoint] = ((clindex&0xff000000)>>24)%18;
3655 xyz[ipoint][0] = GetXrow(row[ipoint]);
3656 xyz[ipoint][1] = cl->GetY();
3657 xyz[ipoint][2] = cl->GetZ();
3661 // Calculate seed state vector and covariance matrix
3663 Double_t alpha, cs,sn, xx2,yy2;
3665 alpha = (sec[1]-sec[2])*fSectors->GetAlpha();
3666 cs = TMath::Cos(alpha);
3667 sn = TMath::Sin(alpha);
3668 xx2= xyz[1][0]*cs-xyz[1][1]*sn;
3669 yy2= xyz[1][0]*sn+xyz[1][1]*cs;
3673 alpha = (sec[0]-sec[2])*fSectors->GetAlpha();
3674 cs = TMath::Cos(alpha);
3675 sn = TMath::Sin(alpha);
3676 xx2= xyz[0][0]*cs-xyz[0][1]*sn;
3677 yy2= xyz[0][0]*sn+xyz[0][1]*cs;
3683 Double_t x[5],c[15];
3687 x[4]=F1(xyz[2][0],xyz[2][1],xyz[1][0],xyz[1][1],xyz[0][0],xyz[0][1]);
3688 x[2]=F2(xyz[2][0],xyz[2][1],xyz[1][0],xyz[1][1],xyz[0][0],xyz[0][1]);
3689 x[3]=F3n(xyz[2][0],xyz[2][1],xyz[0][0],xyz[0][1],xyz[2][2],xyz[0][2],x[4]);
3691 Double_t sy =0.1, sz =0.1;
3693 Double_t sy1=0.2, sz1=0.2;
3694 Double_t sy2=0.2, sz2=0.2;
3697 Double_t f40=(F1(xyz[2][0],xyz[2][1]+sy,xyz[1][0],xyz[1][1],xyz[0][0],xyz[0][1])-x[4])/sy;
3698 Double_t f42=(F1(xyz[2][0],xyz[2][1],xyz[1][0],xyz[1][1]+sy,xyz[0][0],xyz[0][1])-x[4])/sy;
3699 Double_t f43=(F1(xyz[2][0],xyz[2][1],xyz[1][0],xyz[1][1],xyz[0][0],xyz[0][1]+sy)-x[4])/sy;
3700 Double_t f20=(F2(xyz[2][0],xyz[2][1]+sy,xyz[1][0],xyz[1][1],xyz[0][0],xyz[0][1])-x[2])/sy;
3701 Double_t f22=(F2(xyz[2][0],xyz[2][1],xyz[1][0],xyz[1][1]+sy,xyz[0][0],xyz[0][1])-x[2])/sy;
3702 Double_t f23=(F2(xyz[2][0],xyz[2][1],xyz[1][0],xyz[1][1],xyz[0][0],xyz[0][1]+sy)-x[2])/sy;
3704 Double_t f30=(F3(xyz[2][0],xyz[2][1]+sy,xyz[0][0],xyz[0][1],xyz[2][2],xyz[0][2])-x[3])/sy;
3705 Double_t f31=(F3(xyz[2][0],xyz[2][1],xyz[0][0],xyz[0][1],xyz[2][2]+sz,xyz[0][2])-x[3])/sz;
3706 Double_t f32=(F3(xyz[2][0],xyz[2][1],xyz[0][0],xyz[0][1]+sy,xyz[2][2],xyz[0][2])-x[3])/sy;
3707 Double_t f34=(F3(xyz[2][0],xyz[2][1],xyz[0][0],xyz[0][1],xyz[2][2],xyz[0][2]+sz)-x[3])/sz;
3712 c[3]=f20*sy1; c[4]=0.; c[5]=f20*sy1*f20+f22*sy2*f22+f23*sy3*f23;
3713 c[6]=f30*sy1; c[7]=f31*sz1; c[8]=f30*sy1*f20+f32*sy2*f22;
3714 c[9]=f30*sy1*f30+f31*sz1*f31+f32*sy2*f32+f34*sz2*f34;
3715 c[10]=f40*sy1; c[11]=0.; c[12]=f40*sy1*f20+f42*sy2*f22+f43*sy3*f23;
3716 c[13]=f30*sy1*f40+f32*sy2*f42;
3717 c[14]=f40*sy1*f40+f42*sy2*f42+f43*sy3*f43;
3719 // Int_t row1 = fSectors->GetRowNumber(xyz[2][0]);
3720 AliTPCseed *seed=new AliTPCseed(0, x, c, xyz[2][0], sec[2]*fSectors->GetAlpha()+fSectors->GetAlphaShift());
3721 seed->fLastPoint = row[2];
3722 seed->fFirstPoint = row[2];
3726 Int_t AliTPCtrackerMI::CheckKinkPoint(AliTPCseed*seed, Float_t th)
3731 for (Int_t i=0;i<12;i++) seed->fKinkPoint[i]=0;
3733 if (TMath::Abs(seed->GetC())>0.01) return 0;
3736 Float_t x[160], y[160], erry[160], z[160], errz[160];
3738 Float_t xt[160], yt[160], zt[160];
3743 Int_t middle = seed->GetNumberOfClusters()/2;
3746 // find central sector, get local cooordinates
3748 for (Int_t i=seed->fFirstPoint;i<=seed->fLastPoint;i++) {
3749 sec[i]= seed->GetClusterSector(i)%18;
3752 AliTPCclusterMI * cl = seed->fClusterPointer[i];
3753 // if (cl==0) cl = GetClusterMI(seed->GetClusterIndex2(i));
3760 if (i>i2) i2 = i; //last point with cluster
3761 if (i2<i1) i1 = i; //first point with cluster
3764 AliTPCTrackerPoint * point = seed->GetTrackPoint(i);
3766 yt[i] = point->GetY();
3767 zt[i] = point->GetZ();
3769 if (point->GetX()>0){
3770 erry[i] = point->GetErrY();
3771 errz[i] = point->GetErrZ();
3776 secm = sec[i]; //central sector
3777 padm = i; //middle point with cluster
3782 // rotate position to global coordinate system connected to sector at last the point
3784 for (Int_t i=i1;i<=i2;i++){
3786 if (sec[i]<0) continue;
3787 Double_t alpha = (sec[i2]-sec[i])*fSectors->GetAlpha();
3788 Double_t cs = TMath::Cos(alpha);
3789 Double_t sn = TMath::Sin(alpha);
3790 Float_t xx2= x[i]*cs+y[i]*sn;
3791 Float_t yy2= -x[i]*sn+y[i]*cs;
3795 xx2= xt[i]*cs+yt[i]*sn;
3796 yy2= -xt[i]*sn+yt[i]*cs;
3801 //get "state" vector
3802 Double_t xh[5],xm = x[padm];
3805 xh[4]=F1(xt[i2],yt[i2],xt[padm],yt[padm],xt[i1],yt[i1]);
3806 xh[2]=F2(xt[i2],yt[i2],xt[padm],yt[padm],xt[i1],yt[i1]);
3807 xh[3]=F3n(xt[i2],yt[i2],xt[i1],yt[i1],zt[i2],zt[i1],xh[4]);
3810 for (Int_t i=i1;i<=i2;i++){
3812 if (sec[i]<0) continue;
3813 GetProlongation(x[i2], x[i],xh,yy,zz);
3814 if (TMath::Abs(y[i]-yy)>4||TMath::Abs(z[i]-zz)>4){
3816 //xxh[4]=F1old(x[i2],y[i2],x[padm],y[padm],x[i1],y[i1]);
3817 //xxh[2]=F2old(x[i2],y[i2],x[padm],y[padm],x[i1],y[i1]);
3818 Error("AliTPCtrackerMI::CheckKinkPoint","problem\n");
3823 Float_t dyup[160],dydown[160], dzup[160], dzdown[160];
3824 Float_t yup[160], ydown[160], zup[160], zdown[160];
3826 AliTPCpolyTrack ptrack1,ptrack2;
3829 for (Int_t i=i1;i<=i2;i++){
3830 AliTPCclusterMI * cl = seed->fClusterPointer[i];
3832 if (cl->GetType()<0) continue;
3833 if (cl->GetType()>10) continue;
3836 ptrack1.AddPoint(x[i]-xm,y[i],z[i],0.1,0.1);
3838 if (ptrack1.GetN()>4.){
3839 ptrack1.UpdateParameters();
3841 ptrack1.GetFitDerivation(x[i]-xm,ddy,ddz);
3843 ptrack1.GetFitPoint(x[i]-xm,yy,zz);
3852 dyup[i]=0.; //not enough points
3857 for (Int_t i=i2;i>=i1;i--){
3858 AliTPCclusterMI * cl = seed->fClusterPointer[i];
3860 if (cl->GetType()<0) continue;
3861 if (cl->GetType()>10) continue;
3863 ptrack2.AddPoint(x[i]-xm,y[i],z[i],0.1,0.1);
3865 if (ptrack2.GetN()>4){
3866 ptrack2.UpdateParameters();
3868 ptrack2.GetFitDerivation(x[i]-xm,ddy,ddz);
3870 ptrack2.GetFitPoint(x[i]-xm,yy,zz);
3878 dydown[i]=0.; //not enough points
3883 // find maximal difference of the derivation
3884 for (Int_t i=0;i<12;i++) seed->fKinkPoint[i]=0;
3887 for (Int_t i=i1+10;i<i2-10;i++){
3888 if ( (TMath::Abs(dydown[i])<0.00000001) || (TMath::Abs(dyup[i])<0.00000001) ||i<30)continue;
3889 // printf("%f\t%f\t%f\t%f\t%f\n",x[i],dydown[i],dyup[i],dzdown[i],dzup[i]);
3891 Float_t ddy = TMath::Abs(dydown[i]-dyup[i]);
3892 Float_t ddz = TMath::Abs(dzdown[i]-dzup[i]);
3893 if ( (ddy+ddz)> th){
3894 seed->fKinkPoint[0] = i;
3895 seed->fKinkPoint[1] = ddy;
3896 seed->fKinkPoint[2] = ddz;
3903 //write information to the debug tree
3904 TBranch * br = fTreeDebug->GetBranch("debug");
3905 TClonesArray * arr = new TClonesArray("AliTPCTrackPoint2");
3906 arr->ExpandCreateFast(i2-i1);
3907 br->SetAddress(&arr);
3909 AliTPCclusterMI cldummy;
3911 AliTPCTrackPoint2 pdummy;
3912 pdummy.GetTPoint().fIsShared = 10;
3914 Double_t alpha = sec[i2]*fSectors->GetAlpha();
3915 Double_t cs = TMath::Cos(alpha);
3916 Double_t sn = TMath::Sin(alpha);
3918 for (Int_t i=i1;i<i2;i++){
3919 AliTPCTrackPoint2 *trpoint = (AliTPCTrackPoint2*)arr->UncheckedAt(i-i1);
3921 AliTPCclusterMI * cl0 = seed->fClusterPointer[i];
3923 AliTPCTrackerPoint * point = seed->GetTrackPoint(i);
3926 Double_t x = GetXrow(i);
3927 trpoint->GetTPoint() = *point;
3928 trpoint->GetCPoint() = *cl0;
3929 trpoint->GetCPoint().SetQ(TMath::Abs(cl0->GetQ()));
3930 trpoint->fID = seed->GetUniqueID();
3931 trpoint->fLab = seed->GetLabel();
3933 trpoint->fGX = cs *x + sn*point->GetY();
3934 trpoint->fGY = -sn *x + cs*point->GetY() ;
3935 trpoint->fGZ = point->GetZ();
3937 trpoint->fDY = y[i];
3938 trpoint->fDZ = z[i];
3940 trpoint->fDYU = dyup[i];
3941 trpoint->fDZU = dzup[i];
3943 trpoint->fDYD = dydown[i];
3944 trpoint->fDZD = dzdown[i];
3946 if (TMath::Abs(dyup[i])>0.00000000001 &&TMath::Abs(dydown[i])>0.00000000001){
3947 trpoint->fDDY = dydown[i]-dyup[i];
3948 trpoint->fDDZ = dzdown[i]-dzup[i];
3956 trpoint->GetCPoint()= cldummy;
3973 AliTPCseed* AliTPCtrackerMI::ReSeed(AliTPCseed *t)
3976 // reseed - refit - track
3979 // Int_t last = fSectors->GetNRows()-1;
3981 if (fSectors == fOuterSec){
3982 first = TMath::Max(first, t->fFirstPoint-fInnerSec->GetNRows());
3986 first = t->fFirstPoint;
3988 AliTPCseed * seed = MakeSeed(t,0.1,0.5,0.9);
3989 FollowBackProlongation(*t,fSectors->GetNRows()-1);
3991 FollowProlongation(*t,first);
4001 //_____________________________________________________________________________
4002 Int_t AliTPCtrackerMI::ReadSeeds(const TFile *inp) {
4003 //-----------------------------------------------------------------
4004 // This function reades track seeds.
4005 //-----------------------------------------------------------------
4006 TDirectory *savedir=gDirectory;
4008 TFile *in=(TFile*)inp;
4009 if (!in->IsOpen()) {
4010 cerr<<"AliTPCtrackerMI::ReadSeeds(): input file is not open !\n";
4015 TTree *seedTree=(TTree*)in->Get("Seeds");
4017 cerr<<"AliTPCtrackerMI::ReadSeeds(): ";
4018 cerr<<"can't get a tree with track seeds !\n";
4021 AliTPCtrack *seed=new AliTPCtrack;
4022 seedTree->SetBranchAddress("tracks",&seed);
4024 if (fSeeds==0) fSeeds=new TObjArray(15000);
4026 Int_t n=(Int_t)seedTree->GetEntries();
4027 for (Int_t i=0; i<n; i++) {
4028 seedTree->GetEvent(i);
4029 fSeeds->AddLast(new AliTPCseed(*seed,seed->GetAlpha()));
4038 Int_t AliTPCtrackerMI::Clusters2Tracks (AliESD *esd)
4041 if (fSeeds) DeleteSeeds();
4044 if (!fSeeds) return 1;
4051 //_____________________________________________________________________________
4052 Int_t AliTPCtrackerMI::Clusters2Tracks() {
4053 //-----------------------------------------------------------------
4054 // This is a track finder.
4055 //-----------------------------------------------------------------
4056 TDirectory *savedir=gDirectory;
4060 fSeeds = Tracking();
4063 Info("Clusters2Tracks","Time for tracking: \t");timer.Print();timer.Start();
4065 //activate again some tracks
4066 for (Int_t i=0; i<fSeeds->GetEntriesFast(); i++) {
4067 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i), &t=*pt;
4069 Int_t nc=t.GetNumberOfClusters();
4071 delete fSeeds->RemoveAt(i);
4074 if (pt->fRemoval==10) {
4075 if (pt->GetDensityFirst(20)>0.8 || pt->GetDensityFirst(30)>0.8 || pt->GetDensityFirst(40)>0.7)
4076 pt->Desactivate(10); // make track again active
4078 pt->Desactivate(20);
4079 delete fSeeds->RemoveAt(i);
4083 RemoveDouble(fSeeds,0.2,0.6,11);
4084 RemoveUsed(fSeeds,0.5,0.5,6);
4087 Int_t nseed=fSeeds->GetEntriesFast();
4089 for (Int_t i=0; i<nseed; i++) {
4090 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i), &t=*pt;
4092 Int_t nc=t.GetNumberOfClusters();
4094 delete fSeeds->RemoveAt(i);
4097 CookLabel(pt,0.1); //For comparison only
4098 //if ((pt->IsActive() || (pt->fRemoval==10) )&& nc>50 &&pt->GetNumberOfClusters()>0.4*pt->fNFoundable){
4099 if ((pt->IsActive() || (pt->fRemoval==10) )){
4101 if (fDebug>0) cerr<<found<<'\r';
4105 delete fSeeds->RemoveAt(i);
4109 //RemoveOverlap(fSeeds,0.99,7,kTRUE);
4111 //RemoveUsed(fSeeds,0.9,0.9,6);
4113 nseed=fSeeds->GetEntriesFast();
4115 for (Int_t i=0; i<nseed; i++) {
4116 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i), &t=*pt;
4118 Int_t nc=t.GetNumberOfClusters();
4120 delete fSeeds->RemoveAt(i);
4124 t.CookdEdx(0.02,0.6);
4125 // CheckKinkPoint(&t,0.05);
4126 //if ((pt->IsActive() || (pt->fRemoval==10) )&& nc>50 &&pt->GetNumberOfClusters()>0.4*pt->fNFoundable){
4127 if ((pt->IsActive() || (pt->fRemoval==10) )){
4135 delete fSeeds->RemoveAt(i);
4136 //AliTPCseed * seed1 = ReSeed(pt,0.05,0.5,1);
4138 // FollowProlongation(*seed1,0);
4139 // Int_t n = seed1->GetNumberOfClusters();
4140 // printf("fP4\t%f\t%f\n",seed1->GetC(),pt->GetC());
4141 // printf("fN\t%d\t%d\n", seed1->GetNumberOfClusters(),pt->GetNumberOfClusters());
4144 //AliTPCseed * seed2 = ReSeed(pt,0.95,0.5,0.05);
4148 SortTracks(fSeeds, 1);
4152 PrepareForBackProlongation(fSeeds,5.);
4153 PropagateBack(fSeeds);
4154 printf("Time for back propagation: \t");timer.Print();timer.Start();
4158 PrepareForProlongation(fSeeds,5.);
4159 PropagateForward2(fSeeds);
4161 printf("Time for FORWARD propagation: \t");timer.Print();timer.Start();
4162 // RemoveUsed(fSeeds,0.7,0.7,6);
4163 //RemoveOverlap(fSeeds,0.9,7,kTRUE);
4165 nseed=fSeeds->GetEntriesFast();
4167 for (Int_t i=0; i<nseed; i++) {
4168 AliTPCseed *pt=(AliTPCseed*)fSeeds->UncheckedAt(i), &t=*pt;
4170 Int_t nc=t.GetNumberOfClusters();
4172 delete fSeeds->RemoveAt(i);
4175 t.CookdEdx(0.02,0.6);
4176 // CookLabel(pt,0.1); //For comparison only
4177 //if ((pt->IsActive() || (pt->fRemoval==10) )&& nc>50 &&pt->GetNumberOfClusters()>0.4*pt->fNFoundable){
4178 if ((pt->IsActive() || (pt->fRemoval==10) )){
4179 cerr<<found++<<'\r';
4182 delete fSeeds->RemoveAt(i);
4187 // fNTracks = found;
4189 Info("Clusters2Tracks","Time for overlap removal, track writing and dedx cooking: \t"); timer.Print();timer.Start();
4192 // cerr<<"Number of found tracks : "<<"\t"<<found<<endl;
4193 Info("Clusters2Tracks","Number of found tracks %d",found);
4195 // UnloadClusters();
4200 void AliTPCtrackerMI::Tracking(TObjArray * arr)
4203 // tracking of the seeds
4206 fSectors = fOuterSec;
4207 ParallelTracking(arr,150,63);
4208 fSectors = fOuterSec;
4209 ParallelTracking(arr,63,0);
4212 TObjArray * AliTPCtrackerMI::Tracking(Int_t seedtype, Int_t i1, Int_t i2, Float_t cuts[4], Float_t dy, Int_t dsec)
4217 TObjArray * arr = new TObjArray;
4219 fSectors = fOuterSec;
4222 for (Int_t sec=0;sec<fkNOS;sec++){
4223 if (seedtype==3) MakeSeeds3(arr,sec,i1,i2,cuts,dy, dsec);
4224 if (seedtype==4) MakeSeeds5(arr,sec,i1,i2,cuts,dy);
4225 if (seedtype==2) MakeSeeds2(arr,sec,i1,i2,cuts,dy);
4228 Info("Tracking","\nSeeding - %d\t%d\t%d\t%d\n",seedtype,i1,i2,arr->GetEntriesFast());
4240 TObjArray * AliTPCtrackerMI::Tracking()
4246 Int_t nup=fOuterSec->GetNRows()+fInnerSec->GetNRows();
4248 TObjArray * seeds = new TObjArray;
4257 Float_t fnumber = 3.0;
4258 Float_t fdensity = 3.0;
4263 for (Int_t delta = 0; delta<18; delta+=6){
4267 arr = Tracking(3,nup-1-delta,nup-1-delta-gap,cuts,-1,1);
4268 SumTracks(seeds,arr);
4269 SignClusters(seeds,fnumber,fdensity);
4271 for (Int_t i=2;i<6;i+=2){
4272 // seed high pt tracks
4275 arr = Tracking(3,nup-i-delta,nup-i-delta-gap,cuts,-1,0);
4276 SumTracks(seeds,arr);
4277 SignClusters(seeds,fnumber,fdensity);
4282 // RemoveUsed(seeds,0.9,0.9,1);
4283 // UnsignClusters();
4284 // SignClusters(seeds,fnumber,fdensity);
4288 for (Int_t delta = 20; delta<120; delta+=10){
4290 // seed high pt tracks
4294 arr = Tracking(3,nup-delta,nup-delta-gap,cuts,-1);
4295 SumTracks(seeds,arr);
4296 SignClusters(seeds,fnumber,fdensity);
4301 arr = Tracking(3,nup-delta-5,nup-delta-5-gap,cuts,-1);
4302 SumTracks(seeds,arr);
4303 SignClusters(seeds,fnumber,fdensity);
4314 Info("Tracking()","\n\nPrimary seeding\t%d\n\n",seeds->GetEntriesFast());
4318 // RemoveUsed(seeds,0.75,0.75,1);
4320 //SignClusters(seeds,fnumber,fdensity);
4329 arr = Tracking(4,nup-1,nup-1-gap,cuts,-1);
4330 SumTracks(seeds,arr);
4331 SignClusters(seeds,fnumber,fdensity);
4333 arr = Tracking(4,nup-2,nup-2-gap,cuts,-1);
4334 SumTracks(seeds,arr);
4335 SignClusters(seeds,fnumber,fdensity);
4337 arr = Tracking(4,nup-3,nup-3-gap,cuts,-1);
4338 SumTracks(seeds,arr);
4339 SignClusters(seeds,fnumber,fdensity);
4343 for (Int_t delta = 3; delta<30; delta+=5){
4349 arr = Tracking(4,nup-1-delta,nup-1-delta-gap,cuts,-1);
4350 SumTracks(seeds,arr);
4351 SignClusters(seeds,fnumber,fdensity);
4353 arr = Tracking(4,nup-3-delta,nup-5-delta-gap,cuts,4);
4354 SumTracks(seeds,arr);
4355 SignClusters(seeds,fnumber,fdensity);
4367 for (Int_t delta = 30; delta<70; delta+=10){
4373 arr = Tracking(4,nup-1-delta,nup-1-delta-gap,cuts,-1);
4374 SumTracks(seeds,arr);
4375 SignClusters(seeds,fnumber,fdensity);
4377 arr = Tracking(4,nup-5-delta,nup-5-delta-gap,cuts,5 );
4378 SumTracks(seeds,arr);
4379 SignClusters(seeds,fnumber,fdensity);
4383 Info("Tracking()","\n\nSecondary seeding\t%d\n\n",seeds->GetEntriesFast());
4394 void AliTPCtrackerMI::SumTracks(TObjArray *arr1,TObjArray *arr2) const
4397 //sum tracks to common container
4398 //remove suspicious tracks
4399 Int_t nseed = arr2->GetEntriesFast();
4400 for (Int_t i=0;i<nseed;i++){
4401 AliTPCseed *pt=(AliTPCseed*)arr2->UncheckedAt(i);
4404 // NORMAL ACTIVE TRACK
4405 if (pt->IsActive()){
4406 arr1->AddLast(arr2->RemoveAt(i));
4409 //remove not usable tracks
4410 if (pt->fRemoval!=10){
4411 delete arr2->RemoveAt(i);
4414 // REMOVE VERY SHORT TRACKS
4415 if (pt->GetNumberOfClusters()<20){
4416 delete arr2->RemoveAt(i);
4419 // ENABLE ONLY ENOUGH GOOD STOPPED TRACKS
4420 if (pt->GetDensityFirst(20)>0.8 || pt->GetDensityFirst(30)>0.8 || pt->GetDensityFirst(40)>0.7)
4421 arr1->AddLast(arr2->RemoveAt(i));
4423 delete arr2->RemoveAt(i);
4432 void AliTPCtrackerMI::ParallelTracking(TObjArray * arr, Int_t rfirst, Int_t rlast)
4435 // try to track in parralel
4437 Int_t nseed=arr->GetEntriesFast();
4438 //prepare seeds for tracking
4439 for (Int_t i=0; i<nseed; i++) {
4440 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i), &t=*pt;
4442 if (!t.IsActive()) continue;
4443 // follow prolongation to the first layer
4444 if ( (fSectors ==fInnerSec) || (t.fFirstPoint-fParam->GetNRowLow()>rfirst+1) )
4445 FollowProlongation(t, rfirst+1);
4450 for (Int_t nr=rfirst; nr>=rlast; nr--){
4451 if (nr<fInnerSec->GetNRows())
4452 fSectors = fInnerSec;
4454 fSectors = fOuterSec;
4455 // make indexes with the cluster tracks for given
4457 // find nearest cluster
4458 for (Int_t i=0; i<nseed; i++) {
4459 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i), &t=*pt;
4461 if (!pt->IsActive()) continue;
4462 // if ( (fSectors ==fOuterSec) && (pt->fFirstPoint-fParam->GetNRowLow())<nr) continue;
4463 if (pt->fRelativeSector>17) {
4466 UpdateClusters(t,nr);
4468 // prolonagate to the nearest cluster - if founded
4469 for (Int_t i=0; i<nseed; i++) {
4470 AliTPCseed *pt=(AliTPCseed*)arr->UncheckedAt(i);
4472 if (!pt->IsActive()) continue;
4473 // if ((fSectors ==fOuterSec) && (pt->fFirstPoint-fParam->GetNRowLow())<nr) continue;
4474 if (pt->fRelativeSector>17) {
4477 FollowToNextCluster(*pt,nr);
4482 void AliTPCtrackerMI::PrepareForBackProlongation(TObjArray * arr,Float_t fac) const
4486 // if we use TPC track itself we have to "update" covariance
4488 Int_t nseed= arr->GetEntriesFast();
4489 for (Int_t i=0;i<nseed;i++){
4490 AliTPCseed *pt = (AliTPCseed*)arr->UncheckedAt(i);
4494 //rotate to current local system at first accepted point
4495 Int_t index = pt->GetClusterIndex2(pt->fFirstPoint);
4496 Int_t sec = (index&0xff000000)>>24;
4498 Float_t angle1 = fInnerSec->GetAlpha()*sec+fInnerSec->GetAlphaShift();
4499 if (angle1>TMath::Pi())
4500 angle1-=2.*TMath::Pi();
4501 Float_t angle2 = pt->GetAlpha();
4503 if (TMath::Abs(angle1-angle2)>0.001){
4504 pt->Rotate(angle1-angle2);
4505 //angle2 = pt->GetAlpha();
4506 //pt->fRelativeSector = pt->GetAlpha()/fInnerSec->GetAlpha();
4507 //if (pt->GetAlpha()<0)
4508 // pt->fRelativeSector+=18;
4509 //sec = pt->fRelativeSector;
4518 void AliTPCtrackerMI::PrepareForProlongation(TObjArray * arr, Float_t fac) const
4522 // if we use TPC track itself we have to "update" covariance
4524 Int_t nseed= arr->GetEntriesFast();
4525 for (Int_t i=0;i<nseed;i++){
4526 AliTPCseed *pt = (AliTPCseed*)arr->UncheckedAt(i);
4529 pt->fFirstPoint = pt->fLastPoint;
4537 Int_t AliTPCtrackerMI::PropagateBack(TObjArray * arr)
4540 // make back propagation
4542 Int_t nseed= arr->GetEntriesFast();
4543 for (Int_t i=0;i<nseed;i++){
4544 AliTPCseed *pt = (AliTPCseed*)arr->UncheckedAt(i);
4546 //AliTPCseed *pt2 = new AliTPCseed(*pt);
4547 fSectors = fInnerSec;
4548 //FollowBackProlongation(*pt,fInnerSec->GetNRows()-1);
4549 //fSectors = fOuterSec;
4550 FollowBackProlongation(*pt,fInnerSec->GetNRows()+fOuterSec->GetNRows()-1);
4551 //if (pt->GetNumberOfClusters()<(pt->fEsd->GetTPCclusters(0)) ){
4552 // Error("PropagateBack","Not prolonged track %d",pt->GetLabel());
4553 // FollowBackProlongation(*pt2,fInnerSec->GetNRows()+fOuterSec->GetNRows()-1);
4561 Int_t AliTPCtrackerMI::PropagateForward2(TObjArray * arr)
4564 // make forward propagation
4566 Int_t nseed= arr->GetEntriesFast();
4568 for (Int_t i=0;i<nseed;i++){
4569 AliTPCseed *pt = (AliTPCseed*)arr->UncheckedAt(i);
4571 FollowProlongation(*pt,0);
4578 Int_t AliTPCtrackerMI::PropagateForward()
4581 // propagate track forward
4583 Int_t nseed = fSeeds->GetEntriesFast();
4584 for (Int_t i=0;i<nseed;i++){
4585 AliTPCseed *pt = (AliTPCseed*)fSeeds->UncheckedAt(i);
4587 AliTPCseed &t = *pt;
4588 Double_t alpha=t.GetAlpha() - fSectors->GetAlphaShift();
4589 if (alpha > 2.*TMath::Pi()) alpha -= 2.*TMath::Pi();
4590 if (alpha < 0. ) alpha += 2.*TMath::Pi();
4591 t.fRelativeSector = Int_t(alpha/fSectors->GetAlpha()+0.0001)%fN;
4595 fSectors = fOuterSec;
4596 ParallelTracking(fSeeds,fOuterSec->GetNRows()+fInnerSec->GetNRows()-1,fInnerSec->GetNRows());
4597 fSectors = fInnerSec;
4598 ParallelTracking(fSeeds,fInnerSec->GetNRows()-1,0);
4608 Int_t AliTPCtrackerMI::PropagateBack(AliTPCseed * pt, Int_t row0, Int_t row1)
4611 // make back propagation, in between row0 and row1
4615 fSectors = fInnerSec;
4618 if (row1<fSectors->GetNRows())
4621 r1 = fSectors->GetNRows()-1;
4623 if (row0<fSectors->GetNRows()&& r1>0 )
4624 FollowBackProlongation(*pt,r1);
4625 if (row1<=fSectors->GetNRows())
4628 r1 = row1 - fSectors->GetNRows();
4629 if (r1<=0) return 0;
4630 if (r1>=fOuterSec->GetNRows()) return 0;
4631 fSectors = fOuterSec;
4632 return FollowBackProlongation(*pt,r1);
4640 void AliTPCtrackerMI::GetShape(AliTPCseed * seed, Int_t row)
4644 Float_t sd2 = TMath::Abs((fParam->GetZLength()-TMath::Abs(seed->GetZ())))*fParam->GetDiffL()*fParam->GetDiffL();
4645 // Float_t padlength = fParam->GetPadPitchLength(seed->fSector);
4646 Float_t padlength = GetPadPitchLength(row);
4648 Float_t sresy = (seed->fSector < fParam->GetNSector()/2) ? 0.2 :0.3;
4649 Float_t angulary = seed->GetSnp();
4650 angulary = angulary*angulary/(1-angulary*angulary);
4651 seed->fCurrentSigmaY2 = sd2+padlength*padlength*angulary/12.+sresy*sresy;
4653 Float_t sresz = fParam->GetZSigma();
4654 Float_t angularz = seed->GetTgl();
4655 seed->fCurrentSigmaZ2 = sd2+padlength*padlength*angularz*angularz*(1+angulary)/12.+sresz*sresz;
4657 Float_t wy = GetSigmaY(seed);
4658 Float_t wz = GetSigmaZ(seed);
4661 if (TMath::Abs(wy/seed->fCurrentSigmaY2-1)>0.0001 || TMath::Abs(wz/seed->fCurrentSigmaZ2-1)>0.0001 ){
4662 printf("problem\n");
4668 Float_t AliTPCtrackerMI::GetSigmaY(AliTPCseed * seed)
4672 Float_t sd2 = TMath::Abs((fParam->GetZLength()-TMath::Abs(seed->GetZ())))*fParam->GetDiffL()*fParam->GetDiffL();
4673 Float_t padlength = fParam->GetPadPitchLength(seed->fSector);
4674 Float_t sres = (seed->fSector < fParam->GetNSector()/2) ? 0.2 :0.3;
4675 Float_t angular = seed->GetSnp();
4676 angular = angular*angular/(1-angular*angular);
4677 // angular*=angular;
4678 //angular = TMath::Sqrt(angular/(1-angular));
4679 Float_t res = TMath::Sqrt(sd2+padlength*padlength*angular/12.+sres*sres);
4682 Float_t AliTPCtrackerMI::GetSigmaZ(AliTPCseed * seed)
4686 Float_t sd2 = TMath::Abs((fParam->GetZLength()-TMath::Abs(seed->GetZ())))*fParam->GetDiffL()*fParam->GetDiffL();
4687 Float_t padlength = fParam->GetPadPitchLength(seed->fSector);
4688 Float_t sres = fParam->GetZSigma();
4689 Float_t angular = seed->GetTgl();
4690 Float_t res = TMath::Sqrt(sd2+padlength*padlength*angular*angular/12.+sres*sres);
4697 //__________________________________________________________________________
4698 void AliTPCtrackerMI::CookLabel(AliTPCseed *t, Float_t wrong) const {
4699 //--------------------------------------------------------------------
4700 //This function "cooks" a track label. If label<0, this track is fake.
4701 //--------------------------------------------------------------------
4702 Int_t noc=t->GetNumberOfClusters();
4704 //printf("\nnot founded prolongation\n\n\n");
4710 AliTPCclusterMI *clusters[160];
4712 for (Int_t i=0;i<160;i++) {
4719 for (i=0; i<160 && current<noc; i++) {
4721 Int_t index=t->GetClusterIndex2(i);
4722 if (index<=0) continue;
4723 if (index&0x8000) continue;
4725 //clusters[current]=GetClusterMI(index);
4726 if (t->fClusterPointer[i]){
4727 clusters[current]=t->fClusterPointer[i];
4733 Int_t lab=123456789;
4734 for (i=0; i<noc; i++) {
4735 AliTPCclusterMI *c=clusters[i];
4737 lab=TMath::Abs(c->GetLabel(0));
4739 for (j=0; j<noc; j++) if (lb[j]==lab || mx[j]==0) break;
4745 for (i=0; i<noc; i++) if (mx[i]>max) {max=mx[i]; lab=lb[i];}
4747 for (i=0; i<noc; i++) {
4748 AliTPCclusterMI *c=clusters[i];
4750 if (TMath::Abs(c->GetLabel(1)) == lab ||
4751 TMath::Abs(c->GetLabel(2)) == lab ) max++;
4754 if ((1.- Float_t(max)/noc) > wrong) lab=-lab;
4757 Int_t tail=Int_t(0.10*noc);
4760 for (i=1; i<=160&&ind<tail; i++) {
4761 // AliTPCclusterMI *c=clusters[noc-i];
4762 AliTPCclusterMI *c=clusters[i];
4764 if (lab == TMath::Abs(c->GetLabel(0)) ||
4765 lab == TMath::Abs(c->GetLabel(1)) ||
4766 lab == TMath::Abs(c->GetLabel(2))) max++;
4769 if (max < Int_t(0.5*tail)) lab=-lab;
4776 //delete[] clusters;
4780 Int_t AliTPCtrackerMI::AliTPCSector::GetRowNumber(Double_t x) const
4782 //return pad row number for this x
4785 r=fRow[fN-1].GetX();
4786 if (x > r) return fN;
4788 if (x < r) return -1;
4789 return Int_t((x-r)/fPadPitchLength + 0.5);}
4791 r=fRow[fN-1].GetX();
4792 if (x > r) return fN;
4794 if (x < r) return -1;
4795 Double_t r1=fRow[64].GetX();
4797 return Int_t((x-r)/f1PadPitchLength + 0.5);}
4799 return (Int_t((x-r1)/f2PadPitchLength + 0.5)+64);}
4803 //_________________________________________________________________________
4804 void AliTPCtrackerMI::AliTPCSector::Setup(const AliTPCParam *par, Int_t f) {
4805 //-----------------------------------------------------------------------
4806 // Setup inner sector
4807 //-----------------------------------------------------------------------
4809 fAlpha=par->GetInnerAngle();
4810 fAlphaShift=par->GetInnerAngleShift();
4811 fPadPitchWidth=par->GetInnerPadPitchWidth();
4812 fPadPitchLength=par->GetInnerPadPitchLength();
4813 fN=par->GetNRowLow();
4814 fRow=new AliTPCRow[fN];
4815 for (Int_t i=0; i<fN; i++) {
4816 fRow[i].SetX(par->GetPadRowRadiiLow(i));
4817 fRow[i].fDeadZone =1.5; //1.5 cm of dead zone
4820 fAlpha=par->GetOuterAngle();
4821 fAlphaShift=par->GetOuterAngleShift();
4822 fPadPitchWidth = par->GetOuterPadPitchWidth();
4823 fPadPitchLength = par->GetOuter1PadPitchLength();
4824 f1PadPitchLength = par->GetOuter1PadPitchLength();
4825 f2PadPitchLength = par->GetOuter2PadPitchLength();
4827 fN=par->GetNRowUp();
4828 fRow=new AliTPCRow[fN];
4829 for (Int_t i=0; i<fN; i++) {
4830 fRow[i].SetX(par->GetPadRowRadiiUp(i));
4831 fRow[i].fDeadZone =1.5; // 1.5 cm of dead zone
4836 AliTPCtrackerMI::AliTPCRow::AliTPCRow() {
4838 // default constructor
4846 AliTPCtrackerMI::AliTPCRow::~AliTPCRow(){
4853 //_________________________________________________________________________
4855 AliTPCtrackerMI::AliTPCRow::InsertCluster(const AliTPCclusterMI* c, UInt_t index) {
4856 //-----------------------------------------------------------------------
4857 // Insert a cluster into this pad row in accordence with its y-coordinate
4858 //-----------------------------------------------------------------------
4859 if (fN==kMaxClusterPerRow) {
4860 cerr<<"AliTPCRow::InsertCluster(): Too many clusters !\n"; return;
4862 if (fN==0) {fIndex[0]=index; fClusters[fN++]=c; return;}
4863 Int_t i=Find(c->GetZ());
4864 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliTPCclusterMI*));
4865 memmove(fIndex +i+1 ,fIndex +i,(fN-i)*sizeof(UInt_t));
4866 fIndex[i]=index; fClusters[i]=c; fN++;
4869 void AliTPCtrackerMI::AliTPCRow::ResetClusters() {
4875 //delete[] fClusterArray;
4876 if (fClusters1) delete []fClusters1;
4877 if (fClusters2) delete []fClusters2;
4884 //___________________________________________________________________
4885 Int_t AliTPCtrackerMI::AliTPCRow::Find(Double_t z) const {
4886 //-----------------------------------------------------------------------
4887 // Return the index of the nearest cluster
4888 //-----------------------------------------------------------------------
4889 if (fN==0) return 0;
4890 if (z <= fClusters[0]->GetZ()) return 0;
4891 if (z > fClusters[fN-1]->GetZ()) return fN;
4892 Int_t b=0, e=fN-1, m=(b+e)/2;
4893 for (; b<e; m=(b+e)/2) {
4894 if (z > fClusters[m]->GetZ()) b=m+1;
4902 //___________________________________________________________________
4903 AliTPCclusterMI * AliTPCtrackerMI::AliTPCRow::FindNearest(Double_t y, Double_t z, Double_t roady, Double_t roadz) const {
4904 //-----------------------------------------------------------------------
4905 // Return the index of the nearest cluster in z y
4906 //-----------------------------------------------------------------------
4907 Float_t maxdistance = roady*roady + roadz*roadz;
4909 AliTPCclusterMI *cl =0;
4910 for (Int_t i=Find(z-roadz); i<fN; i++) {
4911 AliTPCclusterMI *c=(AliTPCclusterMI*)(fClusters[i]);
4912 if (c->GetZ() > z+roadz) break;
4913 if ( (c->GetY()-y) > roady ) continue;
4914 Float_t distance = (c->GetZ()-z)*(c->GetZ()-z)+(c->GetY()-y)*(c->GetY()-y);
4915 if (maxdistance>distance) {
4916 maxdistance = distance;
4923 AliTPCclusterMI * AliTPCtrackerMI::AliTPCRow::FindNearest2(Double_t y, Double_t z, Double_t roady, Double_t roadz,UInt_t & index) const
4925 //-----------------------------------------------------------------------
4926 // Return the index of the nearest cluster in z y
4927 //-----------------------------------------------------------------------
4928 Float_t maxdistance = roady*roady + roadz*roadz;
4929 Int_t iz1 = TMath::Max(fFastCluster[Int_t(z-roadz+254.5)]-1,0);
4930 Int_t iz2 = TMath::Min(fFastCluster[Int_t(z+roadz+255.5)]+1,fN);
4932 AliTPCclusterMI *cl =0;
4933 //FindNearest3(y,z,roady,roadz,index);
4934 // for (Int_t i=Find(z-roadz); i<fN; i++) {
4935 for (Int_t i=iz1; i<iz2; i++) {
4936 AliTPCclusterMI *c=(AliTPCclusterMI*)(fClusters[i]);
4937 if (c->GetZ() > z+roadz) break;
4938 if ( c->GetY()-y > roady ) continue;
4939 if ( y-c->GetY() > roady ) continue;
4940 Float_t distance = (c->GetZ()-z)*(c->GetZ()-z)+(c->GetY()-y)*(c->GetY()-y);
4941 if (maxdistance>distance) {
4942 maxdistance = distance;
4945 //roady = TMath::Sqrt(maxdistance);
4953 AliTPCclusterMI * AliTPCtrackerMI::AliTPCRow::FindNearest3(Double_t y, Double_t z, Double_t roady, Double_t roadz,UInt_t & index) const
4955 //-----------------------------------------------------------------------
4956 // Return the index of the nearest cluster in z y
4957 //-----------------------------------------------------------------------
4958 Float_t maxdistance = roady*roady + roadz*roadz;
4959 // Int_t iz = Int_t(z+255.);
4960 AliTPCclusterMI *cl =0;
4961 for (Int_t i=Find(z-roadz); i<fN; i++) {
4962 //for (Int_t i=fFastCluster[iz-2]; i<fFastCluster[iz+2]; i++) {
4963 AliTPCclusterMI *c=(AliTPCclusterMI*)(fClusters[i]);
4964 if (c->GetZ() > z+roadz) break;
4965 if ( c->GetY()-y > roady ) continue;
4966 if ( y-c->GetY() > roady ) continue;
4967 Float_t distance = (c->GetZ()-z)*(c->GetZ()-z)+(c->GetY()-y)*(c->GetY()-y);
4968 if (maxdistance>distance) {
4969 maxdistance = distance;
4972 //roady = TMath::Sqrt(maxdistance);
4981 AliTPCseed::AliTPCseed():AliTPCtrack(){
4985 for (Int_t i=0;i<200;i++) SetClusterIndex2(i,-3);
4986 for (Int_t i=0;i<160;i++) fClusterPointer[i]=0;
5003 AliTPCseed::AliTPCseed(const AliTPCseed &s):AliTPCtrack(s){
5004 //---------------------
5005 // dummy copy constructor
5006 //-------------------------
5008 AliTPCseed::AliTPCseed(const AliTPCtrack &t):AliTPCtrack(t){
5017 for (Int_t i=0;i<160;i++) {
5018 fClusterPointer[i] = 0;
5019 Int_t index = t.GetClusterIndex(i);
5021 SetClusterIndex2(i,index);
5024 SetClusterIndex2(i,-3);
5037 AliTPCseed::AliTPCseed(const AliKalmanTrack &t, Double_t a):AliTPCtrack(t,a){
5041 for (Int_t i=0;i<160;i++) {
5042 fClusterPointer[i] = 0;
5043 Int_t index = t.GetClusterIndex(i);
5044 SetClusterIndex2(i,index);
5065 AliTPCseed::AliTPCseed(UInt_t index, const Double_t xx[5], const Double_t cc[15],
5066 Double_t xr, Double_t alpha):
5067 AliTPCtrack(index, xx, cc, xr, alpha) {
5072 for (Int_t i=0;i<200;i++) SetClusterIndex2(i,-3);
5073 for (Int_t i=0;i<160;i++) fClusterPointer[i]=0;
5082 // fHelixIn = new TClonesArray("AliHelix",0);
5083 //fHelixOut = new TClonesArray("AliHelix",0);
5093 AliTPCseed::~AliTPCseed(){
5096 if (fPoints) delete fPoints;
5098 if (fEPoints) delete fEPoints;
5103 AliTPCTrackerPoint * AliTPCseed::GetTrackPoint(Int_t i)
5107 return &fTrackPoints[i];
5110 void AliTPCseed::RebuildSeed()
5113 // rebuild seed to be ready for storing
5114 AliTPCclusterMI cldummy;
5116 AliTPCTrackPoint pdummy;
5117 pdummy.GetTPoint().fIsShared = 10;
5118 for (Int_t i=0;i<160;i++){
5119 AliTPCclusterMI * cl0 = fClusterPointer[i];
5120 AliTPCTrackPoint *trpoint = (AliTPCTrackPoint*)fPoints->UncheckedAt(i);
5122 trpoint->GetTPoint() = *(GetTrackPoint(i));
5123 trpoint->GetCPoint() = *cl0;
5124 trpoint->GetCPoint().SetQ(TMath::Abs(cl0->GetQ()));
5128 trpoint->GetCPoint()= cldummy;
5136 Double_t AliTPCseed::GetDensityFirst(Int_t n)
5140 // return cluster for n rows bellow first point
5141 Int_t nfoundable = 1;
5143 for (Int_t i=fLastPoint-1;i>0&&nfoundable<n; i--){
5144 Int_t index = GetClusterIndex2(i);
5145 if (index!=-1) nfoundable++;
5146 if (index>0) nfound++;
5148 if (nfoundable<n) return 0;
5149 return Double_t(nfound)/Double_t(nfoundable);
5154 void AliTPCseed::GetClusterStatistic(Int_t first, Int_t last, Int_t &found, Int_t &foundable, Int_t &shared, Bool_t plus2)
5156 // get cluster stat. on given region
5161 for (Int_t i=first;i<last; i++){
5162 Int_t index = GetClusterIndex2(i);
5163 if (index!=-1) foundable++;
5164 if (fClusterPointer[i]) {
5170 if (fClusterPointer[i]->IsUsed(10)) {
5174 if (!plus2) continue; //take also neighborhoud
5176 if ( (i>0) && fClusterPointer[i-1]){
5177 if (fClusterPointer[i-1]->IsUsed(10)) {
5182 if ( fClusterPointer[i+1]){
5183 if (fClusterPointer[i+1]->IsUsed(10)) {
5190 //if (shared>found){
5191 //Error("AliTPCseed::GetClusterStatistic","problem\n");
5195 //_____________________________________________________________________________
5196 void AliTPCseed::CookdEdx(Double_t low, Double_t up,Int_t i1, Int_t i2, Bool_t onlyused) {
5197 //-----------------------------------------------------------------
5198 // This funtion calculates dE/dX within the "low" and "up" cuts.
5199 //-----------------------------------------------------------------
5202 Float_t angular[200];
5203 Float_t weight[200];
5206 // TClonesArray & arr = *fPoints;
5207 Float_t meanlog = 100.;
5209 Float_t mean[4] = {0,0,0,0};
5210 Float_t sigma[4] = {1000,1000,1000,1000};
5211 Int_t nc[4] = {0,0,0,0};
5212 Float_t norm[4] = {1000,1000,1000,1000};
5217 for (Int_t of =0; of<4; of++){
5218 for (Int_t i=of+i1;i<i2;i+=4)
5220 Int_t index = fIndex[i];
5221 if (index<0||index&0x8000) continue;
5223 //AliTPCTrackPoint * point = (AliTPCTrackPoint *) arr.At(i);
5224 AliTPCTrackerPoint * point = GetTrackPoint(i);
5225 //AliTPCTrackerPoint * pointm = GetTrackPoint(i-1);
5226 //AliTPCTrackerPoint * pointp = 0;
5227 //if (i<159) pointp = GetTrackPoint(i+1);
5229 if (point==0) continue;
5230 AliTPCclusterMI * cl = fClusterPointer[i];
5231 if (cl==0) continue;
5232 if (onlyused && (!cl->IsUsed(10))) continue;
5233 if (cl->IsUsed(11)) {
5237 Int_t type = cl->GetType();
5238 //if (point->fIsShared){
5243 // if (pointm->fIsShared) continue;
5245 // if (pointp->fIsShared) continue;
5247 if (type<0) continue;
5248 //if (type>10) continue;
5249 //if (point->GetErrY()==0) continue;
5250 //if (point->GetErrZ()==0) continue;
5252 //Float_t ddy = (point->GetY()-cl->GetY())/point->GetErrY();
5253 //Float_t ddz = (point->GetZ()-cl->GetZ())/point->GetErrZ();
5254 //if ((ddy*ddy+ddz*ddz)>10) continue;
5257 // if (point->GetCPoint().GetMax()<5) continue;
5258 if (cl->GetMax()<5) continue;
5259 Float_t angley = point->GetAngleY();
5260 Float_t anglez = point->GetAngleZ();
5262 Float_t rsigmay2 = point->GetSigmaY();
5263 Float_t rsigmaz2 = point->GetSigmaZ();
5267 rsigmay += pointm->GetTPoint().GetSigmaY();
5268 rsigmaz += pointm->GetTPoint().GetSigmaZ();
5272 rsigmay += pointp->GetTPoint().GetSigmaY();
5273 rsigmaz += pointp->GetTPoint().GetSigmaZ();
5280 Float_t rsigma = TMath::Sqrt(rsigmay2*rsigmaz2);
5282 Float_t ampc = 0; // normalization to the number of electrons
5284 // ampc = 1.*point->GetCPoint().GetMax();
5285 ampc = 1.*cl->GetMax();
5286 //ampc = 1.*point->GetCPoint().GetQ();
5287 // AliTPCClusterPoint & p = point->GetCPoint();
5288 // Float_t dy = TMath::Abs(Int_t( TMath::Abs(p.GetY()/0.6)) - TMath::Abs(p.GetY()/0.6)+0.5);
5289 // Float_t iz = (250.0-TMath::Abs(p.GetZ())+0.11)/0.566;
5291 // TMath::Abs( Int_t(iz) - iz + 0.5);
5292 //ampc *= 1.15*(1-0.3*dy);
5293 //ampc *= 1.15*(1-0.3*dz);
5294 // Float_t zfactor = (AliTPCReconstructor::GetCtgRange()-0.0004*TMath::Abs(point->GetCPoint().GetZ()));
5298 //ampc = 1.0*point->GetCPoint().GetMax();
5299 ampc = 1.0*cl->GetMax();
5300 //ampc = 1.0*point->GetCPoint().GetQ();
5301 //AliTPCClusterPoint & p = point->GetCPoint();
5302 // Float_t dy = TMath::Abs(Int_t( TMath::Abs(p.GetY()/0.4)) - TMath::Abs(p.GetY()/0.4)+0.5);
5303 //Float_t iz = (250.0-TMath::Abs(p.GetZ())+0.11)/0.566;
5305 // TMath::Abs( Int_t(iz) - iz + 0.5);
5307 //ampc *= 1.15*(1-0.3*dy);
5308 //ampc *= 1.15*(1-0.3*dz);
5309 // Float_t zfactor = (1.02-0.000*TMath::Abs(point->GetCPoint().GetZ()));
5313 ampc *= 2.0; // put mean value to channel 50
5314 //ampc *= 0.58; // put mean value to channel 50
5316 // if (type>0) w = 1./(type/2.-0.5);
5317 // Float_t z = TMath::Abs(cl->GetZ());
5320 //ampc /= (1+0.0008*z);
5324 //ampc /= (1+0.0008*z);
5326 //ampc /= (1+0.0008*z);
5329 if (type<0) { //amp at the border - lower weight
5334 if (rsigma>1.5) ampc/=1.3; // if big backround
5336 angular[nc[of]] = TMath::Sqrt(1.+angley*angley+anglez*anglez);
5341 TMath::Sort(nc[of],amp,index,kFALSE);
5345 //meanlog = amp[index[Int_t(nc[of]*0.33)]];
5347 for (Int_t i=int(nc[of]*low+0.5);i<int(nc[of]*up+0.5);i++){
5348 Float_t ampl = amp[index[i]]/angular[index[i]];
5349 ampl = meanlog*TMath::Log(1.+ampl/meanlog);
5351 sumw += weight[index[i]];
5352 sumamp += weight[index[i]]*ampl;
5353 sumamp2 += weight[index[i]]*ampl*ampl;
5354 norm[of] += angular[index[i]]*weight[index[i]];
5361 mean[of] = sumamp/sumw;
5362 sigma[of] = sumamp2/sumw-mean[of]*mean[of];
5364 sigma[of] = TMath::Sqrt(sigma[of]);
5368 mean[of] = (TMath::Exp(mean[of]/meanlog)-1)*meanlog;
5369 //mean *=(1-0.02*(sigma/(mean*0.17)-1.));
5370 //mean *=(1-0.1*(norm-1.));
5377 // mean[0]*= (1-0.05*(sigma[0]/(0.01+mean[1]*0.18)-1));
5378 // mean[1]*= (1-0.05*(sigma[1]/(0.01+mean[0]*0.18)-1));
5381 // dedx = (mean[0]* TMath::Sqrt((1.+nc[0]))+ mean[1]* TMath::Sqrt((1.+nc[1])) )/
5382 // ( TMath::Sqrt((1.+nc[0]))+TMath::Sqrt((1.+nc[1])));
5386 for (Int_t i =0;i<4;i++){
5387 if (nc[i]>2&&nc[i]<1000){
5388 dedx += mean[i] *nc[i];
5389 fSdEdx += sigma[i]*(nc[i]-2);
5390 fMAngular += norm[i] *nc[i];
5395 fSDEDX[i] = sigma[i];
5408 // Float_t dedx1 =dedx;
5411 for (Int_t i =0;i<4;i++){
5412 if (nc[i]>2&&nc[i]<1000){
5413 mean[i] = mean[i]*(1-0.12*(sigma[i]/(fSdEdx)-1.));
5414 dedx += mean[i] *nc[i];
5429 Double_t p=TMath::Sqrt((1.+ GetTgl()*GetTgl())/(Get1Pt()*Get1Pt()));
5432 if (dedx < 39.+ 12./(p+0.25)/(p+0.25)) { SetMass(0.13957); return;}
5433 if (dedx < 39.+ 12./p/p) { SetMass(0.49368); return;}
5434 SetMass(0.93827); return;
5438 if (dedx < 39.+ 12./(p+0.25)/(p+0.25)) { SetMass(0.13957); return;}
5439 SetMass(0.93827); return;
5442 SetMass(0.13957); return;
5452 void AliTPCseed::CookdEdx2(Double_t low, Double_t up) {
5453 //-----------------------------------------------------------------
5454 // This funtion calculates dE/dX within the "low" and "up" cuts.
5455 //-----------------------------------------------------------------
5458 Float_t angular[200];
5459 Float_t weight[200];
5461 Bool_t inlimit[200];
5462 for (Int_t i=0;i<200;i++) inlimit[i]=kFALSE;
5463 for (Int_t i=0;i<200;i++) amp[i]=10000;
5464 for (Int_t i=0;i<200;i++) angular[i]= 1;;
5468 Float_t meanlog = 100.;
5469 Int_t indexde[4]={0,64,128,160};
5476 Float_t mean[4] = {0,0,0,0};
5477 Float_t sigma[4] = {1000,1000,1000,1000};
5478 Int_t nc[4] = {0,0,0,0};
5479 Float_t norm[4] = {1000,1000,1000,1000};
5484 // for (Int_t of =0; of<3; of++){
5485 // for (Int_t i=indexde[of];i<indexde[of+1];i++)
5486 for (Int_t i =0; i<160;i++)
5488 AliTPCTrackPoint * point = GetTrackPoint(i);
5489 if (point==0) continue;
5490 if (point->fIsShared){
5494 Int_t type = point->GetCPoint().GetType();
5495 if (type<0) continue;
5496 if (point->GetCPoint().GetMax()<5) continue;
5497 Float_t angley = point->GetTPoint().GetAngleY();
5498 Float_t anglez = point->GetTPoint().GetAngleZ();
5499 Float_t rsigmay = point->GetCPoint().GetSigmaY();
5500 Float_t rsigmaz = point->GetCPoint().GetSigmaZ();
5501 Float_t rsigma = TMath::Sqrt(rsigmay*rsigmaz);
5503 Float_t ampc = 0; // normalization to the number of electrons
5505 ampc = point->GetCPoint().GetMax();
5508 ampc = point->GetCPoint().GetMax();
5510 ampc *= 2.0; // put mean value to channel 50
5511 // ampc *= 0.565; // put mean value to channel 50
5514 Float_t z = TMath::Abs(point->GetCPoint().GetZ());
5521 if (type<0) { //amp at the border - lower weight
5524 if (rsigma>1.5) ampc/=1.3; // if big backround
5525 angular[i] = TMath::Sqrt(1.+angley*angley+anglez*anglez);
5526 amp[i] = ampc/angular[i];
5531 TMath::Sort(159,amp,index,kFALSE);
5532 for (Int_t i=int(anc*low+0.5);i<int(anc*up+0.5);i++){
5533 inlimit[index[i]] = kTRUE; // take all clusters
5536 // meanlog = amp[index[Int_t(anc*0.3)]];
5538 for (Int_t of =0; of<3; of++){
5542 for (Int_t i=indexde[of];i<indexde[of+1];i++)
5544 if (inlimit[i]==kFALSE) continue;
5545 Float_t ampl = amp[i];
5547 ampl = meanlog*TMath::Log(1.+ampl/meanlog);
5550 sumamp += weight[i]*ampl;
5551 sumamp2 += weight[i]*ampl*ampl;
5552 norm[of] += angular[i]*weight[i];
5560 mean[of] = sumamp/sumw;
5561 sigma[of] = sumamp2/sumw-mean[of]*mean[of];
5563 sigma[of] = TMath::Sqrt(sigma[of]);
5566 mean[of] = (TMath::Exp(mean[of]/meanlog)-1)*meanlog;
5576 Float_t www[3] = {12.,14.,17.};
5577 //Float_t www[3] = {1.,1.,1.};
5579 for (Int_t i =0;i<3;i++){
5580 if (nc[i]>2&&nc[i]<1000){
5581 dedx += mean[i] *nc[i]*www[i]/sigma[i];
5582 fSdEdx += sigma[i]*(nc[i]-2)*www[i]/sigma[i];
5583 fMAngular += norm[i] *nc[i];
5584 norm2 += nc[i]*www[i]/sigma[i];
5585 norm3 += (nc[i]-2)*www[i]/sigma[i];
5588 fSDEDX[i] = sigma[i];
5601 // Float_t dedx1 =dedx;
5605 for (Int_t i =0;i<3;i++){
5606 if (nc[i]>2&&nc[i]<1000&&sigma[i]>3){
5607 //mean[i] = mean[i]*(1+0.08*(sigma[i]/(fSdEdx)-1.));
5608 dedx += mean[i] *(nc[i])/(sigma[i]);
5609 norm4 += (nc[i])/(sigma[i]);
5613 if (norm4>0) dedx /= norm4;